Composed Map#

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

Provides a definition of a transport map defined by the composition of several maps.

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

\[ T(x) = T_L(x) \circ \ \dots \ \circ T_1(x), \]
where each map \(T_i:\mathbb{R}^{N}\rightarrow \mathbb{R}^{N}\) is a square transport map.

Checkpointing can be used for deep maps to reduce the amount of memory that is required for gradient computations. The checkpointing logic was adapted from the revolve algorithm of [Griewank and Walther, 2000], which is an optimal binomial checkpointing scheme.

Public Functions

ComposedMap(std::vector<std::shared_ptr<ConditionalMapBase<MemorySpace>>> const &maps, bool moveCoeffs = false, int maxChecks = -1)#

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

Parameters:
  • maps – A vector of ConditionalMapBase objects defining each \(T_k\) in the composition. Note: each map must be square (i.e., have equal input and output dimension).

  • moveCoeffs – Whether you want this ComposedMap object to be responsible for the coefficients already in maps. If true, the new object will take ownership of all coefficient vectors within maps. If false, the new object will not have any coefficients set.

  • maxChecks – The maximum number of checkpoints to use during gradient computations. If maxChecks==1, then no checkpointing will be utilized and all forward states will be recomputed. If maxChecks==components.size(), then all states will be stored and reused during the backward pass. This is the most efficient option, but can require an intractable amount of memory for high-dimensional or deep parameterizations. The default value is -1, which will set the maximum number of checkpoints to be equal to the number of layers (i.e., map.size()).

virtual ~ComposedMap() = 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 InverseImpl(StridedMatrix<const double, MemorySpace> const &x1, StridedMatrix<const double, MemorySpace> const &r, StridedMatrix<double, MemorySpace> output) override#

Pure abstract function overridden by child classes.

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#
virtual void GradientImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) override#