Parameterized Function Base#

template<typename MemorySpace>
class ParameterizedFunctionBase : public std::enable_shared_from_this<ParameterizedFunctionBase<MemorySpace>>#

Subclassed by mpart::AffineFunction< MemorySpace >, mpart::ConditionalMapBase< MemorySpace >, mpart::MultivariateExpansion< BasisEvaluatorType, MemorySpace >

Public Functions

inline ParameterizedFunctionBase(unsigned int inDim, unsigned int outDim, unsigned int nCoeffs)#

Construct a new parameterized function object.

Parameters:
  • inDim – The dimension \(N\) of the input to this map.

  • outDim – The dimension \(M\) of the output from this map.

  • nCoeffs – The number of coefficients in the parameterization.

inline ParameterizedFunctionBase(unsigned int inDim, unsigned int outDim, unsigned int nCoeffs, Kokkos::View<const double*, MemorySpace> coeffsIn)#
virtual ~ParameterizedFunctionBase() = default#
inline virtual Kokkos::View<double*, MemorySpace> &Coeffs()#

Returns a view containing the coefficients for this conditional map. This function returns a reference and can therefore be used to to update the coefficients or even set them to be a subview of a larger view. When the values of the larger view are updated, the subview stored by this class will also be updated. This is particularly convenient when simultaneously optimizing the coefficients over many conditional maps because each map can just use a slice into the larger vector of all coefficients that is updated by the optimizer.

Returns:

A reference to the Kokkos::View containing the coefficients.

virtual void SetCoeffs(Kokkos::View<const double*, MemorySpace> coeffs)#

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.

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

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.

virtual void SetCoeffs(Eigen::Ref<Eigen::VectorXd> coeffs)#

SetCoeffs function with conversion from Eigen to Kokkos vector types.

void WrapCoeffs(Eigen::Ref<Eigen::VectorXd> coeffs)#
virtual Eigen::Map<Eigen::VectorXd> CoeffMap()#

Returns an eigen map wrapping around the coefficient vector, which is stored in a Kokkos::View. Updating the components of this map should also update the view.

inline virtual Kokkos::View<const double*, MemorySpace> Coeffs() const#

Const version of the Coeffs() function.

template<typename ViewType>
inline StridedMatrix<double, typename ViewType::memory_space> Evaluate(ViewType pts)#

Evaluate function with conversion between default view layout and const strided matrix.

Eigen::RowMatrixXd Evaluate(Eigen::Ref<const Eigen::RowMatrixXd> const &pts)#

Evaluate function with conversion from Eigen to Kokkos (and possibly copy to/from device.)

template<typename AnyMemorySpace>
StridedMatrix<double, AnyMemorySpace> Evaluate(StridedMatrix<const double, AnyMemorySpace> const &pts)#

Evaluate the function at multiple points.

Computes \(\mathbf{y}^{(i)}=T(\mathbf{x}^{(i)}; \mathbf{w})\) for \(N\) points \(\{\mathbf{x}^{(1)},\ldots,\mathbf{x}^{(N)}\}\). The parameters \(\mathbf{w}\) are defined by the :code:SetCoeffs function.

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.

Returns:

A \(d_{out}\times N\) matrix containing evaluations of this function.

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

Pure virtual EvaluateImpl function that is overridden by children of this class.

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.

  • output – A preallocated \(d_{out}\times N\) that will be filled with the evaluations. This matrix must be sized correctly before being passed into this function.

template<typename ViewType1, typename ViewType2>
inline StridedMatrix<double, typename ViewType1::memory_space> Gradient(ViewType1 pts, ViewType2 sens)#

Evaluate the gradient of the function with conversion between default view layout and const strided matrix.

Eigen::RowMatrixXd Gradient(Eigen::Ref<const Eigen::RowMatrixXd> const &pts, Eigen::Ref<const Eigen::RowMatrixXd> const &sens)#

Evaluate the gradient of the function with conversion from Eigen to Kokkos (and possibly copy to/from device.)

template<typename AnyMemorySpace>
StridedMatrix<double, AnyMemorySpace> Gradient(StridedMatrix<const double, AnyMemorySpace> const &pts, StridedMatrix<const double, AnyMemorySpace> const &sens)#

Evaluate the gradient of the function at multiple points.

For input points \(x^{(i)}\) and sensitivity vectors \(s^{(i)}\), this function computes

\[ g^{(i)} = \left[s^{(i)}\right]^T\nabla_x T(x^{(i)}; w), \]
where \(\nabla_x T\) is the Jacobian of the function \(T\) with respect to the input \(x\). Note that this function can be used to evaluate one step of the chain rule (e.g., one backpropagation step). Given a scalar-valued function \(g\), the gradient of \(g(T(x))\) with respect to \(x\) is given by \(\left(\nabla g\right)^T \left(\nabla_x T\right)\). Passing \(\nabla g\) as the sensitivity input to this function then allows you to compute this product and thus the gradient of the composed function \(g(T(x))\).

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

  • sens – A \(d_{out}\times N\) matrix containing \(N\) sensitivity vectors in \(\mathbb{R}^{d_{out}}\).

Returns:

A \(d_{in}\times N\) matrix containing the gradient of vectors. Each column corresponds to the gradient at a particular point and sensitivity.

virtual void GradientImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) = 0#
template<typename AnyMemorySpace>
StridedMatrix<double, AnyMemorySpace> CoeffGrad(StridedMatrix<const double, AnyMemorySpace> const &pts, StridedMatrix<const double, AnyMemorySpace> const &sens)#

Computes the gradient of the map output with respect to the map coefficients.

Consider a map \(T(x; w) : \mathbb{R}^N \rightarrow \mathbb{R}^M\) parameterized by coefficients \(w\in\mathbb{R}^K\). This function computes

\[ g_i = s_i^T\nabla_w T(x_i; w), \]
where \(\nabla_w T(x_i; w)\in\mathbb{R}^{M\times K}\) denotes the Jacobian of the map \(T(x_i; w)\) with respect to the coefficients at the single point \(x_i\) and \(s_i\in \mathbb{R}^M\) is a vector of sensitivities. Often the sensititivities \(s_i\) represent the gradient of some objective function with respect to the map output, i.e. \(s_i = \nabla_{y_i} c(y_i)\), where \(c:\mathbb{R}^M\rightarrow \mathbb{R}\) is a scalar-valued objective function and \(y_i=T(x_i;w)\) is the output of the map. In this setting, the vector \(g_i\) computed by this function represents \(g_i = \nabla_{w} c(T(x_i; w))\); and this function essentially computes a single step in the chain rule.

Parameters:
  • pts – A collection of points \(x_i\) where we want to compute the Jacobian. Each column contains a single point.

  • sens – A collection of sensitivity vectors \(s_i\) for each point. Each column is a single \(s_i\) vector and this view should therefore have the same number of columns as pts. It should also have \(M\) rows.

Returns:

A collection of vectors \(g_i\). Will have the same number of columns as pts with \(K\) rows.

template<typename PtsViewType, typename SensViewType>
inline StridedMatrix<double, typename PtsViewType::memory_space> CoeffGrad(PtsViewType pts, SensViewType sens)#

CoeffGrad function with conversion between general view type and const strided matrix.

Eigen::RowMatrixXd CoeffGrad(Eigen::Ref<const Eigen::RowMatrixXd> const &pts, Eigen::Ref<const Eigen::RowMatrixXd> const &sens)#

Coeff grad function with conversion between Eigen and Kokkos matrix types.

virtual void CoeffGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) = 0#

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.

bool CheckCoefficients() const#

Checks to see if the coefficients have been initialized yet, returns true if so, false if not

Public Members

const unsigned int inputDim#
const unsigned int outputDim#

The total dimension of the input N+M.

const unsigned int numCoeffs#

The output dimension M.