Monotone Component#

template<class ExpansionType, class PosFuncType, class QuadratureType, typename MemorySpace>
class MonotoneComponent : public mpart::ConditionalMapBase<MemorySpace>#

Defines a function \(T:R^N\rightarrow R\) such that \(\partial T / \partial x_N >0\) is strictly positive.

The function \(T\) is based on another (generally non-monotone) function \(f : R^N\rightarrow R\), a strictly positve function \(g : R\rightarrow R_{>0}\), and a non-negative scalar \(\epsilon\). Together, these functions define the monotone component $T$ through

\[ T(x_1, x_2, ..., x_D) = f(x_1,x_2,..., x_{D-1}, 0) + \int_0^{x_D} g\left( \partial_D f(x_1,x_2,..., x_{D-1}, t) + \epsilon \right) dt \]

Template Parameters:
  • ExpansionType – A class defining the function \(f\). It must satisfy the cached parameterization concept.

  • PosFuncType – A class defining the function \(g\). This class must have Evaluate and Derivative functions accepting a double and returning a double. The MParT::SoftPlus and MParT::Exp classes in PositiveBijectors.h are examples of classes defining this interface.

  • QuadratureType – A class defining the integration scheme used to approximate \(\int_0^{x_N} g\left( \frac{\partial f}{\partial x_d}(x_1,x_2,..., x_{N-1}, t) \right) dt\). The type must have a function Integrate(f,lb,ub) that accepts a functor f, a double lower bound lb, a double upper bound ub, and returns a double with an estimate of the integral. The MParT::AdaptiveSimpson and MParT::RecursiveQuadrature classes provide this interface.

Public Functions

inline MonotoneComponent(ExpansionType const &expansion, QuadratureType const &quad, bool useContDeriv = true, double nugget = 0.0)#

Construct a monotone component with a specific expansion and quadrature type.

Parameters:
  • expansion – The expansion used to define the function \(f\).

  • quad – The quadrature rule used to approximate \(\int_0^{x_D} g\left( \partial_D f(x_1,x_2,..., x_{D-1}, t) \right) dt\)

  • useContDeriv – A flag to specify whether the analytic derivative of \(T(x_1, x_2, ..., x_D)\) should be used by default, or if the derivative of the discretized integral should be used. If “true”, the analytic or “continuous” derivative will be used. If “false”, the derivative of the numerically approximated integral will be used.

  • nugget – The value of the non negative quantity \(\epsilon\) defining the minimum slope of the transport map.

    See the Diagonal Derivatives mathematical background for more details.

inline MonotoneComponent(ExpansionType const &expansion, QuadratureType const &quad, bool useContDeriv, double nugget, Kokkos::View<const double*, MemorySpace> coeffsIn)#
inline virtual std::shared_ptr<ParameterizedFunctionBase<MemorySpace>> GetBaseFunction() override#

For Monotone parameterizations that are based on a non-monotone base function, this function will return the base function. If the monotone parameterization is not constructed from a non-monotone base, then this function will return a nullptr.

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

Override the ConditionalMapBase Evaluate function.

inline 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.

inline virtual void LogDeterminantImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<double, MemorySpace> output) override#
inline bool isGradFunctionInputValid(int sensRows, int sensCols, int ptsRows, int ptsCols, int outputRows, int outputCols, int expectedOutputRows)#
inline void checkGradFunctionInput(const std::string method, int sensRows, int sensCols, int ptsRows, int ptsCols, int outputRows, int outputCols, int expectedOutputRows)#
inline virtual void GradientImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<const double, MemorySpace> const &sens, StridedMatrix<double, MemorySpace> output) override#
inline 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.

inline virtual void LogDeterminantCoeffGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<double, MemorySpace> output) override#
inline virtual void LogDeterminantInputGradImpl(StridedMatrix<const double, MemorySpace> const &pts, StridedMatrix<double, MemorySpace> output) override#
template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space, class ...OtherTraits>
inline void EvaluateImpl(Kokkos::View<double**, OtherTraits...> const &pts, StridedVector<double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> output)#

Support calling EvaluateImpl with non-const views.

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space, class ...OtherTraits>
inline void EvaluateImpl(Kokkos::View<const double**, OtherTraits...> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> output)#

Evaluates the monotone function \(T(x_1,\ldots,x_D)\) at multiple points.

Parameters:
  • pts[in] A \(D\times N\) array containing the \(N\) points in \(\mathbb{R}^D\) where we want to evaluate the monotone component. Each column is a point.

  • coeffs[in] The coefficients in the expansion defining \(f\). The length of this array must be the same as the number of terms in the multiindex set passed to the constructor.

  • output[out] Kokkos::View<double*> An array containing the evaluattions \(T(x^{(i)}_1,\ldots,x^{(i)}_D)\) for each \(i\in\{0,\ldots,N\}\).

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void InverseImpl(StridedMatrix<const double, MemorySpace> const &xs, StridedVector<const double, MemorySpace> const &ys, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> output, std::map<std::string, std::string> options = std::map<std::string, std::string>())#

Evaluates the inverse of the diagonal of the monotone component.

This function solve the nonlinear equation \(y_D = T(x_1,\ldots,x_D)\) for \(x_D\) given \(N\) different values of \(y_D^{(1)},\ldots,y_D^{(N)}\). It is possible to specify either a single \(x_{1:D-1}\) sample or \(N\) different values \(x_{1:D-1}^{(1)},\ldots,x_{1:D-1}^{(N)}\). If only a single point is given, then this function returns a vector containing \(x_D^{(i)}=T^{-1}(x_1,\ldots,x_{D-1}, y^{(i)}_D)\). When \(N\)N points are given, this function returns a vector containing \(x_D^{(i)}=T^{-1}(x_1^{(i)},\ldots,x_{D-1}^{(i)}, y^{(i)}_D)\).

Parameters:
  • xs – A \(D\times N_1\) array containing \(N_1\) \(x_{1:D-1}\) points. Note that this matrix must have either 1 or N columns and at least \(D\) rows. The \(D\) row will serve as an initial guess for \(x_D\) during the inversion.

  • ys – A length \(N\) array containing \(N\) values of \(y_D\) for use in the solve.

  • coeffs – The coefficients in the expansion defining \(f\). The length of this array must be the same as the number of terms in the multiindex set passed to the constructor.

  • output – An array for storing the computed values of \(y_D^{(i)}\). Memory for this array must be preallocated before calling this function.

  • options – A map containing options for the method (e.g., converge criteria, step sizes). Available options are “Method” (must be “Bracket”), “xtol” (any nonnegative float), and “ytol” (any nonnegative float).

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline Kokkos::View<double*, MemorySpace> ContinuousDerivative(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs)#

Approximates the “continuous derivative” \(\frac{\partial T}{\partial x_D}\) derived from the exact integral form of the transport map.

See the mathematical background section for more details on discrete and continuous map derivatives.

Parameters:
  • pts – The points where we want to evaluate the derivative. Each column represents one point.

  • coeffs – The ceofficients in an expansion for \(f\).

Returns:

Kokkos::View<double*> The values of \(\frac{\partial T}{\partial x_D}\) at each point.

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void ContinuousDerivative(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> derivs)#

Approximates the “continuous derivative” \(\frac{\partial T}{\partial x_D}\) derived from the exact integral form of the transport map.

See the mathematical background section for more details on discrete and continuous map derivatives.

Parameters:
  • pts[in] The points where we want to evaluate the derivative. Each column represents one point.

  • coeffs[in] The ceofficients in an expansion for \(f\).

  • evals[inout] The values of map component itself \(T\) at each point.

  • derivs[inout] The values of \( \frac{\partial T}{\partial x_D}\) at each point.

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline Kokkos::View<double*, MemorySpace> DiscreteDerivative(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs)#

Approximates the “discrete derivative” of the quadrature-based approximation \(\tilde{T}\).

See the mathematical background section for more details on discrete and continuous map derivatives.

Parameters:
  • pts[in] The points where we want to evaluate the derivative. Each column represents one point.

  • coeffs[in] The ceofficients in an expansion for \(f\).

Returns:

Kokkos::View<double*> The values of \( \frac{\partial \tilde{T}}{\partial x_D}\) at each point.

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void DiscreteDerivative(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> evals, StridedVector<double, MemorySpace> derivs)#

Approximates the “discrete derivative” of the quadrature-based approximation \(\tilde{T}\).

See the mathematical background section for more details on discrete and continuous map derivatives.

Parameters:
  • pts[in] The points where we want to evaluate the derivative. Each column represents one point.

  • coeffs[in] The ceofficients in an expansion for \(f\).

  • evals[out] The values of the map component itself \(\tilde{T}\) at each point.

  • derivs[out] Kokkos::View<double*> The values of \( \frac{\partial \tilde{T}}{\partial x_D}\) at each point.

inline bool isJacobianInputValid(int jacRows, int jacCols, int evalRows, int expectJacRows, int expectJacCols, int expectEvalRows)#
inline void checkJacobianInput(const std::string method, int jacRows, int jacCols, int evalRows, int expectJacRows, int expectJacCols, int expectEvalRows)#
template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void CoeffJacobian(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> evaluations, StridedMatrix<double, MemorySpace> jacobian)#

Returns the gradient of the map with respect to the parameters \(\mathbf{w}\) at multiple points.

Consider \(N\) points \(\{\mathbf{x}^{(1)},\ldots,\mathbf{x}^{(N)}\}\) and let \(y_d^{(i)} = T_d(\mathbf{x}^{(i)}; \mathbf{w})\). This function computes \(\nabla_{\mathbf{w}} y_d^{(i)}\) for each output \(y_d^{(i)}\).

See also

CoeffGradient

Parameters:
  • pts[in] A \(D\times N\) matrix containing the points \(x^{(1)},\ldots,x^{(N)}\). Each column is a point.

  • coeffs[in] A vector of coefficients defining the function \(f(\mathbf{x}; \mathbf{w})\).

  • evaluations[out] A vector containing the \(N\) predictions \(y_d^{(i)}\). The vector must be preallocated and have \(N\) components when passed to this function. An error will occur if this vector is not the correct size.

  • jacobian[out] A matrix containing the \(M\times N\) Jacobian matrix, where \(M\) is the length of the parameter vector \(\mathbf{w}\). This matrix must be sized correctly or an error will occur.

template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void InputJacobian(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedVector<double, MemorySpace> evaluations, StridedMatrix<double, MemorySpace> jacobian)#

Returns the gradient of the map with respect to the input \(x_{1:d}\) at multiple points.

Consider \(N\) points \(\{\mathbf{x}^{(1)},\ldots,\mathbf{x}^{(N)}\}\) and let \(y_d^{(i)} = T_d(\mathbf{x}^{(i)}; \mathbf{w})\). This function computes \(\nabla_{\mathbf{x}} y_d^{(i)}\) for each output \(y_d^{(i)}\).

See also

CoeffGradient

Parameters:
  • pts[in] A \(D\times N\) matrix containing the points \(x^{(1)},\ldots,x^{(N)}\). Each column is a point.

  • coeffs[in] A vector of coefficients defining the function \(f(\mathbf{x}; \mathbf{w})\).

  • evaluations[out] A vector containing the \(N\) predictions \(y_d^{(i)}\). The vector must be preallocated and have \(N\) components when passed to this function. An error will occur if this vector is not the correct size.

  • jacobian[out] A matrix containing the \(d\times N\) Jacobian matrix, where \(d\) is the length of the input vector \(\mathbf{x}\). This matrix must be sized correctly or an error will occur.

inline void checkMixedJacobianInput(const std::string method, int jacRows, int jacCols, int expectJacRows, int expectJacCols)#
template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void ContinuousMixedJacobian(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedMatrix<double, MemorySpace> jacobian)#
template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void ContinuousMixedInputJacobian(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedMatrix<double, MemorySpace> jacobian)#
template<typename ExecutionSpace = typename MemoryToExecution<MemorySpace>::Space>
inline void DiscreteMixedJacobian(StridedMatrix<const double, MemorySpace> const &pts, StridedVector<const double, MemorySpace> const &coeffs, StridedMatrix<double, MemorySpace> jacobian)#
inline FixedMultiIndexSet<MemorySpace> GetMultiIndexSet() const#

Give access to the underlying FixedMultiIndexSet

Returns:

The FixedMultiIndexSet

inline std::vector<unsigned int> DiagonalCoeffIndices() const#

Indices of coeffs that correspond to terms non-constant in \(x_d\).

Returns:

std::vector<unsigned int>

Public Static Functions

template<typename PointType, typename CoeffsType>
static inline double EvaluateSingle(double *cache, double *workspace, PointType const &pt, double xd, CoeffsType const &coeffs, QuadratureType const &quad, ExpansionType const &expansion, double nugget = 0.0)#

Evaluates the monotone component at a single point using an existing cache.

Template Parameters:

PointType – The type of point used (subview, etc…)

Parameters:
  • cache – Memory used by the MultivariateExpansion to cache evaluations

  • workspace – Memory used by the quadrature routine to store evaluations

  • pt

  • coeffs

Returns:

double