MUQ  0.4.3
muq::Approximation::SmolyakEstimator< EstimateType > Class Template Referenceabstract

#include <SmolyakEstimator.h>

Detailed Description

template<typename EstimateType>
class muq::Approximation::SmolyakEstimator< EstimateType >

Definition at line 16 of file SmolyakEstimator.h.

Classes

struct  SmolyTerm
 

Public Member Functions

 SmolyakEstimator (std::shared_ptr< muq::Modeling::ModPiece > const &modelIn)
 
virtual ~SmolyakEstimator ()=default
 
virtual EstimateType Compute (std::shared_ptr< muq::Utilities::MultiIndexSet > const &fixedSet, boost::property_tree::ptree options=boost::property_tree::ptree())
 
virtual EstimateType Adapt (boost::property_tree::ptree options)
 
virtual double Error () const
 
virtual unsigned int NumEvals () const
 
virtual std::vector< double > ErrorHistory () const
 
virtual std::vector< int > EvalHistory () const
 
virtual std::vector< double > TimeHistory () const
 
virtual std::vector< std::vector< Eigen::VectorXd > > PointHistory () const
 
virtual std::vector< std::vector< std::shared_ptr< muq::Utilities::MultiIndex > > > TermHistory () const
 

Constructor & Destructor Documentation

◆ SmolyakEstimator()

template<typename EstimateType >
SmolyakEstimator::SmolyakEstimator ( std::shared_ptr< muq::Modeling::ModPiece > const &  modelIn)

Definition at line 13 of file SmolyakEstimator.cpp.

◆ ~SmolyakEstimator()

template<typename EstimateType >
virtual muq::Approximation::SmolyakEstimator< EstimateType >::~SmolyakEstimator ( )
virtualdefault

Member Function Documentation

◆ Adapt()

template<typename EstimateType >
EstimateType SmolyakEstimator::Adapt ( boost::property_tree::ptree  options)
virtual

To be called after Compute, this will continue to refine the estimate until the stopping criteria in options is met.

Definition at line 73 of file SmolyakEstimator.cpp.

◆ AddEstimates()

template<typename EstimateType >
virtual EstimateType muq::Approximation::SmolyakEstimator< EstimateType >::AddEstimates ( double  w1,
EstimateType const &  part1,
double  w2,
EstimateType const &  part2 
) const
protectedpure virtual

◆ AddTerms() [1/2]

template<typename EstimateType >
void SmolyakEstimator::AddTerms ( std::shared_ptr< muq::Utilities::MultiIndexSet > const &  fixedSet)
protectedvirtual

Definition at line 163 of file SmolyakEstimator.cpp.

◆ AddTerms() [2/2]

template<typename EstimateType >
virtual void muq::Approximation::SmolyakEstimator< EstimateType >::AddTerms ( std::vector< std::shared_ptr< muq::Utilities::MultiIndex >> const &  fixedSet)
protectedvirtual

◆ AddToCache()

template<typename EstimateType >
int SmolyakEstimator::AddToCache ( Eigen::VectorXd const &  newPt)
protected

Definition at line 389 of file SmolyakEstimator.cpp.

◆ CacheSize()

template<typename EstimateType >
int muq::Approximation::SmolyakEstimator< EstimateType >::CacheSize ( ) const
inlineprotected

◆ Compute()

template<typename EstimateType >
EstimateType SmolyakEstimator::Compute ( std::shared_ptr< muq::Utilities::MultiIndexSet > const &  fixedSet,
boost::property_tree::ptree  options = boost::property_tree::ptree() 
)
virtual

This is the main function to constructing static or adaptive Smolyak estimates.

Definition at line 324 of file SmolyakEstimator.cpp.

◆ ComputeMagnitude()

template<typename EstimateType >
virtual double muq::Approximation::SmolyakEstimator< EstimateType >::ComputeMagnitude ( EstimateType const &  estimate) const
protectedpure virtual

◆ ComputeOneTerm()

template<typename EstimateType >
virtual EstimateType muq::Approximation::SmolyakEstimator< EstimateType >::ComputeOneTerm ( std::shared_ptr< muq::Utilities::MultiIndex > const &  multi,
std::vector< std::reference_wrapper< const Eigen::VectorXd >> const &  modEvals 
)
protectedpure virtual

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.

Implemented in muq::Approximation::AdaptiveSmolyakQuadrature, and muq::Approximation::AdaptiveSmolyakPCE.

◆ ComputeWeightedSum() [1/2]

template<typename EstimateType >
EstimateType SmolyakEstimator::ComputeWeightedSum
protectedvirtual

Definition at line 50 of file SmolyakEstimator.cpp.

◆ ComputeWeightedSum() [2/2]

template<typename EstimateType >
EstimateType SmolyakEstimator::ComputeWeightedSum ( Eigen::VectorXd const &  weights) const
protectedvirtual

Should compute sum(smolyVals[i] * smolyWeights[i]) and return the result

Definition at line 24 of file SmolyakEstimator.cpp.

◆ Error()

template<typename EstimateType >
virtual double muq::Approximation::SmolyakEstimator< EstimateType >::Error ( ) const
inlinevirtual

Returns the current estimate of the global error in the Smolyak approximation.

Definition at line 38 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::globalError.

◆ ErrorHistory()

template<typename EstimateType >
virtual std::vector<double> muq::Approximation::SmolyakEstimator< EstimateType >::ErrorHistory ( ) const
inlinevirtual

Returns the history of the global error for each adaptation iteration.

Definition at line 44 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::errorHistory.

◆ EvalHistory()

template<typename EstimateType >
virtual std::vector<int> muq::Approximation::SmolyakEstimator< EstimateType >::EvalHistory ( ) const
inlinevirtual

Returns the cumulative number of evaluations after each adaptation iteration. Note that when Compute() is called, the number of evaluations is reset but the evaluation cache is not. If Compute() is called more than once, this may lead to counterintuitive evaluation histories because the points that were evaluated in the first call to Compute() will not be added to the history.

Definition at line 52 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::evalHistory.

◆ EvaluatePoints()

template<typename EstimateType >
void SmolyakEstimator::EvaluatePoints ( std::set< unsigned int > const &  ptsToEval)
protectedvirtual

Evaluates the model at specified points in the cache and saves the results to the evalCache vector.

Parameters
[in]ptsToEvalA set of indices into the cache with points that need evaluation.

Definition at line 299 of file SmolyakEstimator.cpp.

◆ GetFromCache()

template<typename EstimateType >
Eigen::VectorXd const& muq::Approximation::SmolyakEstimator< EstimateType >::GetFromCache ( unsigned int  index) const
inlineprotected

◆ InCache()

template<typename EstimateType >
int muq::Approximation::SmolyakEstimator< EstimateType >::InCache ( Eigen::VectorXd const &  input) const
protected

◆ NumEvals()

template<typename EstimateType >
virtual unsigned int muq::Approximation::SmolyakEstimator< EstimateType >::NumEvals ( ) const
inlinevirtual

Returns the number of model evaluations this factory has performed.

Definition at line 41 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::numEvals.

◆ OneTermPoints()

template<typename EstimateType >
virtual std::vector<Eigen::VectorXd> muq::Approximation::SmolyakEstimator< EstimateType >::OneTermPoints ( std::shared_ptr< muq::Utilities::MultiIndex > const &  multi)
protectedpure virtual

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.

Implemented in muq::Approximation::AdaptiveSmolyakQuadrature, and muq::Approximation::AdaptiveSmolyakPCE.

◆ PointHistory()

template<typename EstimateType >
std::vector< std::vector< Eigen::VectorXd > > SmolyakEstimator::PointHistory
virtual

Return the points that were evaluated during each adaptation step.

Definition at line 173 of file SmolyakEstimator.cpp.

◆ Refine()

template<typename EstimateType >
bool SmolyakEstimator::Refine
protectedvirtual

Refine the approximation by adding to the computed terms. Returns true if any new terms were added.

Definition at line 111 of file SmolyakEstimator.cpp.

◆ Reset()

template<typename EstimateType >
void SmolyakEstimator::Reset
protectedvirtual

Definition at line 358 of file SmolyakEstimator.cpp.

◆ TermHistory()

template<typename EstimateType >
virtual std::vector<std::vector<std::shared_ptr<muq::Utilities::MultiIndex> > > muq::Approximation::SmolyakEstimator< EstimateType >::TermHistory ( ) const
inlinevirtual

Returns the terms used in the estimate as the a function of adaptation iterations.

Definition at line 61 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::termHistory.

◆ TimeHistory()

template<typename EstimateType >
virtual std::vector<double> muq::Approximation::SmolyakEstimator< EstimateType >::TimeHistory ( ) const
inlinevirtual

Returns the cumulative runtime (in seconds) for each adaptation iteration.

Definition at line 55 of file SmolyakEstimator.h.

References muq::Approximation::SmolyakEstimator< EstimateType >::timeHistory.

◆ UpdateErrors()

template<typename EstimateType >
void SmolyakEstimator::UpdateErrors
protectedvirtual

Updates the local error indicators for terms on the leading edge.

Definition at line 308 of file SmolyakEstimator.cpp.

Member Data Documentation

◆ cacheTol

template<typename EstimateType >
const double muq::Approximation::SmolyakEstimator< EstimateType >::cacheTol = 10.0*std::numeric_limits<double>::epsilon()
protected

Definition at line 146 of file SmolyakEstimator.h.

◆ errorHistory

template<typename EstimateType >
std::vector<double> muq::Approximation::SmolyakEstimator< EstimateType >::errorHistory
protected

Holds the history of the error. Each component corresponds to an iteration.

Definition at line 123 of file SmolyakEstimator.h.

Referenced by muq::Approximation::SmolyakEstimator< EstimateType >::ErrorHistory().

◆ errorTol

template<typename EstimateType >
double muq::Approximation::SmolyakEstimator< EstimateType >::errorTol = 0.0
protected

Tolerance on the global error indicator to continue adapting.

Definition at line 152 of file SmolyakEstimator.h.

◆ evalCache

template<typename EstimateType >
std::vector<Eigen::VectorXd> muq::Approximation::SmolyakEstimator< EstimateType >::evalCache
protected

Definition at line 139 of file SmolyakEstimator.h.

◆ evalHistory

template<typename EstimateType >
std::vector<int> muq::Approximation::SmolyakEstimator< EstimateType >::evalHistory
protected

Holds the history of how many model evaluations have occured. Each component corresponds to an adaptation iteration.

Definition at line 126 of file SmolyakEstimator.h.

Referenced by muq::Approximation::SmolyakEstimator< EstimateType >::EvalHistory().

◆ globalError

template<typename EstimateType >
double muq::Approximation::SmolyakEstimator< EstimateType >::globalError
protected

◆ maxNumEvals

template<typename EstimateType >
unsigned int muq::Approximation::SmolyakEstimator< EstimateType >::maxNumEvals = std::numeric_limits<unsigned int>::max()
protected

Tolerance on the maximum number of evaluations.

Definition at line 155 of file SmolyakEstimator.h.

◆ model

template<typename EstimateType >
std::shared_ptr<muq::Modeling::ModPiece> muq::Approximation::SmolyakEstimator< EstimateType >::model
protected

The model used to construct the approximations.

Definition at line 117 of file SmolyakEstimator.h.

◆ numEvals

template<typename EstimateType >
unsigned int muq::Approximation::SmolyakEstimator< EstimateType >::numEvals =0
protected

The number of model evaluations that have been performed.

Definition at line 158 of file SmolyakEstimator.h.

Referenced by muq::Approximation::SmolyakEstimator< EstimateType >::NumEvals().

◆ nzTol

template<typename EstimateType >
const double muq::Approximation::SmolyakEstimator< EstimateType >::nzTol = 10.0*std::numeric_limits<double>::epsilon()
protected

Definition at line 197 of file SmolyakEstimator.h.

◆ pointCache

template<typename EstimateType >
muq::Modeling::DynamicKDTreeAdaptor muq::Approximation::SmolyakEstimator< EstimateType >::pointCache
protected

◆ pointHistory

template<typename EstimateType >
std::vector<std::set<unsigned int> > muq::Approximation::SmolyakEstimator< EstimateType >::pointHistory
protected

Indices in the cache for the points that were evaluated during each adaptation iteration.

Definition at line 132 of file SmolyakEstimator.h.

◆ termHistory

template<typename EstimateType >
std::vector<std::vector<std::shared_ptr<muq::Utilities::MultiIndex> > > muq::Approximation::SmolyakEstimator< EstimateType >::termHistory
protected

The terms that are added during each refinement.

Definition at line 135 of file SmolyakEstimator.h.

Referenced by muq::Approximation::SmolyakEstimator< EstimateType >::TermHistory().

◆ termMultis

template<typename EstimateType >
std::shared_ptr<muq::Utilities::MultiIndexSet> muq::Approximation::SmolyakEstimator< EstimateType >::termMultis
protected

Multiindices defining each tensor product term in the Smolyak approximation.

Definition at line 120 of file SmolyakEstimator.h.

◆ terms

template<typename EstimateType >
std::vector<SmolyTerm> muq::Approximation::SmolyakEstimator< EstimateType >::terms
protected

Definition at line 193 of file SmolyakEstimator.h.

◆ timeHistory

template<typename EstimateType >
std::vector<double> muq::Approximation::SmolyakEstimator< EstimateType >::timeHistory
protected

Holds the history of how. Each component corresponds to an adaptation iteration.

Definition at line 129 of file SmolyakEstimator.h.

Referenced by muq::Approximation::SmolyakEstimator< EstimateType >::TimeHistory().

◆ timeTol

template<typename EstimateType >
double muq::Approximation::SmolyakEstimator< EstimateType >::timeTol = std::numeric_limits<double>::infinity()
protected

Tolerance on the time (in seconds) allowed to continue adapting.

Definition at line 149 of file SmolyakEstimator.h.


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