Triangular Map#

template<typename MemorySpace>
class TriangularMap : public mpart::ConditionalMapBase<MemorySpace>#

Provides a definition of block lower triangular transport maps.

This class defines a map \(T:\mathbb{R}^N\rightarrow \mathbb{R}^M\) with the block triangular structure

\[\begin{split} T(x) = \left[\begin{array}{l} T_1(x_{1:N_1})\\ \vdots \\ T_k(x_{1:N_2})\\ \vdots T_K(x_{1:N}) \end{array} \right], \end{split}\]
where each component \(T_i(x_{1:N_i}):\mathbb{R}^{N_i}\rightarrow \mathbb{R}^{M_i}\) is a function depending on the first \(N_i\) inputs and returning \(M_i\) outputs. Note that this function must be invertible in the last \(M_i\) input arguments. For example, it must be possible to solve \(T_i(x_1:{N_i-M_i}, x_{N_i-M_i:N_i}) = r\) for \(x_{N_i-M_i:N_i}\) given a vector \(r\in\mathbb{R}^{M_i}\) and previous components \(x_1:{N_i-M_i}\).

This block triangular form is analogous to a block triangular matrix where each \(M_i\timesM_i\) diagonal block is positive definite.

Public Functions

TriangularMap(std::vector<std::shared_ptr<ConditionalMapBase<MemorySpace>>> const &components, bool moveCoeffs = false)#

Construct a block triangular map from a collection of other ConditionalMapBase objects.

Parameters:
  • components – A vector of ConditionalMapBase objects defining each \(T_k\) in the block triangular map. To maintain the correct block structure, the dimensions of the components must satisfy \(N_k = N_{k-1}+M_{k}\).

  • moveCoeffs – A boolean to choose whether to keep the coefficients from the maps in components or make new ones. If true, the new object takes ownership of all the coefficients vectors within all the maps in components (changing the coefficients in the new map will then change the coefficients in the original maps). If false, no coefficients are copied or created.

virtual ~TriangularMap() = default#
virtual void SetCoeffs(Kokkos::View<const double*, MemorySpace> coeffs) override#

Set the internally stored view of coefficients.

Performs a deep copy of the input coefficients to the internally stored coefficients.

Parameters:

coeffs – A view containing the coefficients to copy.

virtual void WrapCoeffs(Kokkos::View<double*, MemorySpace> coeffs) override#

Wrap the internal coefficient view around another view.

Performs a shallow copy of the input coefficients to the internally stored coefficients. If values in the view passed to this function are changed, the values will also change in the internally stored view.

Parameters:

coeffs – A view containing the coefficients we want to wrap.

inline virtual std::shared_ptr<ConditionalMapBase<MemorySpace>> GetComponent(unsigned int i)#
virtual void LogDeterminantImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<double, MemorySpace> output) override#

Computes the log determinant of the Jacobian matrix of this map.

Parameters:
  • pts – The points where the jacobian should be evaluated. Should have \(N\) rows. Each column contains a single point.

  • output – A vector with length equal to the number of columns in pts containing the log determinant of the Jacobian. This vector should be correctly allocated and sized before calling this function.

virtual void EvaluateImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<double, MemorySpace> output) override#

Evaluates the map.

Parameters:
  • pts – The points where the map should be evaluated. Should have \(N\) rows. Each column contains a single point.

  • output – A matrix with \(M\) rows to store the map output. Should have the same number of columns as pts and be allocated before calling this function.

virtual void GradientImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) override#
virtual void InverseImpl(StridedMatrix<const double, MemorySpace> const &x1, StridedMatrix<const double, MemorySpace> const &r, StridedMatrix<double, MemorySpace> output) override#

Evaluates the map inverse.

To understand this function, consider splitting the map input \(x_{1:N}\) into two parts so that \(x_{1:N} = [x_{1:N-M},x_{N-M+1:M}]\). Note that the second part of this split \(x_{N-M+1:M}\) has the same dimension as the map output. With this split, the map becomes \(T(x_{1:N-M},x_{N-M+1:M})=r_{1:M}\). Given \(x_{1:N-M}\) and \(r_{1:M}\), the InverseImpl function solves for the value of \(x_{N-M+1:M}\) that satisfies \(T(x_{1:N-M},x_{N-M+1:M})=r_{1:M}\). Our shorthand for this will be \(x_{N-M+1:M} = T^{-1}(r_{1:M}; x_{1:N-M})\).

Parameters:
  • x1 – The points \(x_{1:N-M}\) where the map inverse should be evaluated. Each column contains a single point.

  • r – The map output \(r_{1:M}\) to invert.

  • output – A matrix with \(M\) rows to store the computed map inverse \(x_{N-M+1:M}\).

virtual void InverseInplace(StridedMatrix<double, MemorySpace> x1, StridedMatrix<const double, MemorySpace> const &r)#
virtual void CoeffGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) override#

Pure virtual function overridden by child classes for computing the gradient of the function output with respect to the parameter vector \(\mathbf{w}\). See the non-virtual CoeffGrad function for more details.

Evaluates the gradient with respect to the coefficients and stores the results in a preallocated matrix.

Parameters:
  • pts – A \(d_{in}\times N\) matrix containining \(N\) points in \(\mathbb{R}^d\) where this function be evaluated. Each column is a point.

  • sens – A matrix of sensitivity vectors. Each column contains one sensitivity.

  • output – A preallocated matrix to store the results.

virtual void LogDeterminantCoeffGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<double, MemorySpace> output) override#
virtual void LogDeterminantInputGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<double, MemorySpace> output) override#
std::vector<unsigned int> DiagonalCoeffIndices() const#