MUQ  0.4.3
muq::Approximation::AdaptiveSmolyakQuadrature Class Reference

#include <AdaptiveSmolyakQuadrature.h>

Inheritance diagram for muq::Approximation::AdaptiveSmolyakQuadrature:

Detailed Description

Definition at line 12 of file AdaptiveSmolyakQuadrature.h.

Public Member Functions

 AdaptiveSmolyakQuadrature (std::shared_ptr< muq::Modeling::ModPiece > const &modelIn, std::vector< std::shared_ptr< Quadrature >> const &quad1d)
 
virtual ~AdaptiveSmolyakQuadrature ()=default
 

Constructor & Destructor Documentation

◆ AdaptiveSmolyakQuadrature()

AdaptiveSmolyakQuadrature::AdaptiveSmolyakQuadrature ( std::shared_ptr< muq::Modeling::ModPiece > const &  modelIn,
std::vector< std::shared_ptr< Quadrature >> const &  quad1d 
)

Definition at line 8 of file AdaptiveSmolyakQuadrature.cpp.

◆ ~AdaptiveSmolyakQuadrature()

virtual muq::Approximation::AdaptiveSmolyakQuadrature::~AdaptiveSmolyakQuadrature ( )
virtualdefault

Member Function Documentation

◆ AddEstimates()

Eigen::VectorXd AdaptiveSmolyakQuadrature::AddEstimates ( double  w1,
Eigen::VectorXd const &  part1,
double  w2,
Eigen::VectorXd const &  part2 
) const
overrideprotectedvirtual

◆ ComputeMagnitude()

double AdaptiveSmolyakQuadrature::ComputeMagnitude ( Eigen::VectorXd const &  estimate) const
overrideprotectedvirtual

◆ ComputeOneTerm()

Eigen::VectorXd AdaptiveSmolyakQuadrature::ComputeOneTerm ( std::shared_ptr< muq::Utilities::MultiIndex > const &  multi,
std::vector< std::reference_wrapper< const Eigen::VectorXd >> const &  modEvals 
)
overrideprotectedvirtual

This function works in tandem with the OneTermPoints function. After the model has been evaluated at the points returned by OneTermPoints, this function will compute an estimate with the new model evaluations. In the quadrature setting, the estimate will be computed with a weighted sum of the model evaluations stored in the modEvals vector.

Implements muq::Approximation::SmolyakEstimator< Eigen::VectorXd >.

Definition at line 34 of file AdaptiveSmolyakQuadrature.cpp.

References cachedMulti, muq::Approximation::FullTensorQuadrature::Compute(), tensQuad, and muq::Approximation::Quadrature::Weights().

◆ OneTermPoints()

std::vector< Eigen::VectorXd > AdaptiveSmolyakQuadrature::OneTermPoints ( std::shared_ptr< muq::Utilities::MultiIndex > const &  multi)
overrideprotectedvirtual

Computes the locations where the model will need to be evaluated in order to construct a single tensor-product estimate. For example, in the Smolyak quadrature setting, this function will return the quadrature points coming from the tensor product quadrature rule defined by the multiindex.

This function works in tandem with the ComputeOneTerm function, which takes model evaluations and actually returns the tensor product estimate. These two functions are split to allow the SmolyakEstimator to handle any job scheduling or caching that may be needed for parallel model evaluations.

Implements muq::Approximation::SmolyakEstimator< Eigen::VectorXd >.

Definition at line 21 of file AdaptiveSmolyakQuadrature.cpp.

References cachedMulti, muq::Approximation::FullTensorQuadrature::Compute(), muq::Approximation::Quadrature::Points(), and tensQuad.

Member Data Documentation

◆ cachedMulti

std::shared_ptr<muq::Utilities::MultiIndex> muq::Approximation::AdaptiveSmolyakQuadrature::cachedMulti
protected

Definition at line 32 of file AdaptiveSmolyakQuadrature.h.

Referenced by ComputeOneTerm(), and OneTermPoints().

◆ tensQuad

FullTensorQuadrature muq::Approximation::AdaptiveSmolyakQuadrature::tensQuad
protected

Definition at line 33 of file AdaptiveSmolyakQuadrature.h.

Referenced by ComputeOneTerm(), and OneTermPoints().


The documentation for this class was generated from the following files: