SG++-Doxygen-Documentation
|
Interface class for simple usage of the combigrid module. More...
#include <CombigridMultiOperation.hpp>
Public Member Functions | |
CombigridMultiOperation (std::vector< std::shared_ptr< AbstractPointHierarchy >> pointHierarchies, std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> evaluatorPrototypes, std::shared_ptr< LevelManager > levelManager, MultiFunction func, bool exploitNesting=true, FullGridSummationStrategyType summationStrategyType=FullGridSummationStrategyType::LINEAR, std::shared_ptr< NormStrategy< FloatArrayVector >> normStrategy=nullptr) | |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and function to evaluate. More... | |
CombigridMultiOperation (std::vector< std::shared_ptr< AbstractPointHierarchy >> pointHierarchies, std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> evaluatorPrototypes, std::shared_ptr< LevelManager > levelManager, std::shared_ptr< AbstractCombigridStorage > storage, FullGridSummationStrategyType summationStrategyType=FullGridSummationStrategyType::LINEAR, std::shared_ptr< NormStrategy< FloatArrayVector >> normStrategy=nullptr) | |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and CombigridStorage that contains the function to evaluate. More... | |
CombigridMultiOperation (std::vector< std::shared_ptr< AbstractPointHierarchy >> pointHierarchies, std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> evaluatorPrototypes, std::shared_ptr< LevelManager > levelManager, GridFunction gridFunc, bool exploitNesting=true, FullGridSummationStrategyType summationStrategyType=FullGridSummationStrategyType::LINEAR, std::shared_ptr< NormStrategy< FloatArrayVector >> normStrategy=nullptr) | |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and grid function together with the information whether the same point on different levels should be able to have different function values. More... | |
base::DataVector | evaluate (size_t q, std::vector< base::DataVector > const ¶ms) |
Computes the result with regular levels up to 1-norm q (levels start from zero) and with parameter params. More... | |
base::DataVector | evaluate (size_t q, base::DataMatrix const ¶ms=base::DataMatrix(0, 0)) |
See the other version of evaluate() and the setParameters() overloads. More... | |
std::shared_ptr< AbstractMultiStorage< FloatArrayVector > > | getDifferences () |
Returns a storage of the differences (Deltas) that have been computed by the CombigridEvaluator. More... | |
std::shared_ptr< LevelManager > | getLevelManager () |
Via the LevelManager, more options are available than are provided directly by this class. More... | |
std::vector< std::shared_ptr< AbstractPointHierarchy > > | getPointHierarchies () |
base::DataVector | getResult () |
Returns the current computed value. More... | |
std::shared_ptr< AbstractCombigridStorage > | getStorage () |
size_t | getUpperPointBound () const |
size_t | numDims () |
size_t | numGridPoints () |
size_t | numStoredFunctionValues () |
void | setLevelManager (std::shared_ptr< LevelManager > levelManager) |
Can be used to set the level manager, e.g. More... | |
void | setParameters (std::vector< base::DataVector > const ¶ms) |
Sets the parameters for upcoming computations and clears the data structures (removes old computed data). More... | |
void | setParameters (base::DataMatrix const ¶ms) |
Sets the parameters for upcoming computations and clears the data structures (removes old computed data). More... | |
Static Public Member Functions | |
static std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > | createBsplineLinearCoefficientOperation (size_t degree, size_t numDimensions, std::shared_ptr< sgpp::combigrid::AbstractCombigridStorage > coefficientStorage) |
creates a B spline interpolation operation from a storage of interpolation coefficients More... | |
static std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > | createBsplineLinearRefinementOperation (size_t degree, size_t numDimensions, sgpp::combigrid::MultiFunction func, std::shared_ptr< sgpp::combigrid::LevelManager > levelManager) |
Creates a level structure according to an averaging level manager using linear calculations on each level as norm. More... | |
static std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > | createBsplineVarianceRefinementOperation (size_t degree, size_t numDimensions, sgpp::combigrid::MultiFunction func, std::shared_ptr< sgpp::combigrid::LevelManager > levelManager, sgpp::combigrid::WeightFunctionsCollection weightFunctions, sgpp::base::DataVector bounds) |
Creates a level structure according to an averaging level manager using variance calculations on each level as norm. More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpChebyshevPolynomialInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing polynomial interpolation on a Chebyshev grid with an exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpClenshawCurtisPolynomialInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with an exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpClenshawCurtisQuadrature (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a Clenshaw-Curtis with exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpL2LejaPolynomialInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with an exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpLejaPolynomialInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with an exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpUniformBoundaryBsplineInterpolation (size_t numDimensions, MultiFunction func, size_t degree=3) |
Returns a CombigridMultiOperation doing B-spline interpolation on a uniform grid with boundary using an exponential growth strategy (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpUniformBoundaryBsplineQuadrature (size_t numDimensions, MultiFunction func, size_t degree) |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Bsplines) on a uniform grid with exponential growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpUniformBoundaryBsplineSquareQuadrature (size_t numDimensions, MultiFunction func, size_t degree) |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Bsplines) of f^2, where f is the objective function. More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpUniformBoundaryLinearInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid with boundary using an exponential growth strategy (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createExpUniformLinearInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid without boundary using an exponential growth strategy (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createLinearClenshawCurtisPolynomialInterpolation (size_t numDimensions, MultiFunction func) |
Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with a linear growth (not nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createLinearL2LejaPolynomialInterpolation (size_t numDimensions, MultiFunction func, size_t growthFactor=2) |
Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with linear growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createLinearL2LejaQuadrature (size_t numDimensions, MultiFunction func, size_t growthFactor=2) |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a L2Leja grid with linear growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createLinearLejaPolynomialInterpolation (size_t numDimensions, MultiFunction func, size_t growthFactor=2) |
Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with linear growth (nested points). More... | |
static std::shared_ptr< CombigridMultiOperation > | createLinearLejaQuadrature (size_t numDimensions, MultiFunction func, size_t growthFactor=2) |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a Leja grid with linear growth (nested points). More... | |
Interface class for simple usage of the combigrid module.
Via a CombigridMultiOperation, the evaluation (at multiple interpolation points at once) of the computation pipeline can be easily managed. There are two main ways to create this class:
Via the LevelManager, which can be get and set, one can control which adaptivity criterion might be used. For easy evaluation, there is an evaluate()-method, which does all the work at once and generates a regular level structure. To get more control over the level structure, one may proceed as follows:
sgpp::combigrid::CombigridMultiOperation::CombigridMultiOperation | ( | std::vector< std::shared_ptr< AbstractPointHierarchy >> | pointHierarchies, |
std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> | evaluatorPrototypes, | ||
std::shared_ptr< LevelManager > | levelManager, | ||
MultiFunction | func, | ||
bool | exploitNesting = true , |
||
FullGridSummationStrategyType | summationStrategyType = FullGridSummationStrategyType::LINEAR , |
||
std::shared_ptr< NormStrategy< FloatArrayVector >> | normStrategy = nullptr |
||
) |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and function to evaluate.
References sgpp::combigrid::CombigridMultiOperationImpl::levelManager, and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
sgpp::combigrid::CombigridMultiOperation::CombigridMultiOperation | ( | std::vector< std::shared_ptr< AbstractPointHierarchy >> | pointHierarchies, |
std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> | evaluatorPrototypes, | ||
std::shared_ptr< LevelManager > | levelManager, | ||
std::shared_ptr< AbstractCombigridStorage > | storage, | ||
FullGridSummationStrategyType | summationStrategyType = FullGridSummationStrategyType::LINEAR , |
||
std::shared_ptr< NormStrategy< FloatArrayVector >> | normStrategy = nullptr |
||
) |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and CombigridStorage that contains the function to evaluate.
References sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies, and sgpp::combigrid::CombigridMultiOperationImpl::storage.
sgpp::combigrid::CombigridMultiOperation::CombigridMultiOperation | ( | std::vector< std::shared_ptr< AbstractPointHierarchy >> | pointHierarchies, |
std::vector< std::shared_ptr< AbstractLinearEvaluator< FloatArrayVector >>> | evaluatorPrototypes, | ||
std::shared_ptr< LevelManager > | levelManager, | ||
GridFunction | gridFunc, | ||
bool | exploitNesting = true , |
||
FullGridSummationStrategyType | summationStrategyType = FullGridSummationStrategyType::LINEAR , |
||
std::shared_ptr< NormStrategy< FloatArrayVector >> | normStrategy = nullptr |
||
) |
Constructs a CombigridMultiOperation with the given hierarchies, evaluators, level manager and grid function together with the information whether the same point on different levels should be able to have different function values.
References sgpp::combigrid::CombigridMultiOperationImpl::levelManager, and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
|
static |
creates a B spline interpolation operation from a storage of interpolation coefficients
degree | degree of the b-splines |
numDimensions | number of dimensions |
coefficientStorage | storage of the b-spline coefficients |
References sgpp::combigrid::CombiEvaluators::createCombiMultiEvaluator(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::LINEAR, sgpp::combigrid::Multi_BSplineInterpolation, and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
Referenced by sgpp::combigrid::BsplineStochasticCollocation::updateConfig().
|
static |
Creates a level structure according to an averaging level manager using linear calculations on each level as norm.
This is a very specific case created for the CO2 example. There it serves as a dummy for storing the levelstructure
degree | B spline degree |
numDimensions | number of dimensions |
func | the objective function |
levelManager | level manager |
References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiEvaluators::createCombiMultiEvaluator(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::LINEAR, sgpp::combigrid::Multi_BSplineInterpolation, and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
|
static |
Creates a level structure according to an averaging level manager using variance calculations on each level as norm.
This is a very specific case created for the CO2 example. It can (should?) be generalized
degree | B spline degree |
numDimensions | number of dimensions |
func | the objective function |
levelManager | level manager |
weightFunctions | weight functions |
bounds | bounding box |
References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiEvaluators::createCombiMultiEvaluator(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::Multi_BSplineScalarProduct, sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies, and sgpp::combigrid::VARIANCE.
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a Chebyshev grid with an exponential growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expChebyshev(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with an exponential growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expClenshawCurtis(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a Clenshaw-Curtis with exponential growth (nested points).
Note: This method is not useful as a MultiOperation because the quadrature does not need any parameters. Use CombigridOperation instead.
numDimensions | Dimensionality of the problem. |
func | Function to be integrated. |
References sgpp::combigrid::CombiHierarchies::expClenshawCurtis(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with an exponential growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expL2Leja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with an exponential growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expLeja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing B-spline interpolation on a uniform grid with boundary using an exponential growth strategy (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
degree | degree of the B-spline basis functions |
References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombiEvaluators::multiBSplineInterpolation(), and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
|
static |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Bsplines) on a uniform grid with exponential growth (nested points).
Note: This method is not useful as a MultiOperation because the quadrature does not need any parameters. Use CombigridOperation instead.
numDimensions | Dimensionality of the problem. |
func | Function to be integrated. |
degree | B-spline degree |
References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombiEvaluators::multiBSplineQuadrature(), and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.
|
static |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Bsplines) of f^2, where f is the objective function.
This is needed for variance calculations. on a uniform grid with exponential growth (nested points). Note: This method is not useful as a MultiOperation because the quadrature does not need any parameters. Use CombigridOperation instead.
numDimensions | Dimensionality of the problem. |
func | Function to be integrated. |
degree | B-spline degree |
References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiEvaluators::BSplineScalarProduct(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies, and sgpp::combigrid::QUADRATIC.
|
static |
Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid with boundary using an exponential growth strategy (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expUniformBoundary(), and sgpp::combigrid::CombiEvaluators::multiLinearInterpolation().
|
static |
Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid without boundary using an exponential growth strategy (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiHierarchies::expUniform(), and sgpp::combigrid::CombiEvaluators::multiLinearInterpolation().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with a linear growth (not nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
References sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with linear growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
growthFactor | Parameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used. |
References sgpp::combigrid::CombiHierarchies::linearL2Leja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a L2Leja grid with linear growth (nested points).
Note: This method is not useful as a MultiOperation because the quadrature does not need any parameters. Use CombigridOperation instead.
numDimensions | Dimensionality of the problem. |
func | Function to be integrated. |
growthFactor | Parameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used. |
References sgpp::combigrid::CombiHierarchies::linearL2Leja(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().
|
static |
Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with linear growth (nested points).
numDimensions | Dimensionality of the problem. |
func | Function to be interpolated. |
growthFactor | Parameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used. |
References sgpp::combigrid::CombiHierarchies::linearLeja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().
|
static |
Returns a CombigridMultiOperation doing quadrature (based on integrals of Lagrange polynomials) on a Leja grid with linear growth (nested points).
Note: This method is not useful as a MultiOperation because the quadrature does not need any parameters. Use CombigridOperation instead.
numDimensions | Dimensionality of the problem. |
func | Function to be integrated. |
growthFactor | Parameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used. |
References sgpp::combigrid::CombiHierarchies::linearLeja(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().
base::DataVector sgpp::combigrid::CombigridMultiOperation::evaluate | ( | size_t | q, |
std::vector< base::DataVector > const & | params | ||
) |
Computes the result with regular levels up to 1-norm q (levels start from zero) and with parameter params.
This is a convenience function. If you need other functionality, use getLevelManager() and operate directly on the LevelManager.
base::DataVector sgpp::combigrid::CombigridMultiOperation::evaluate | ( | size_t | q, |
base::DataMatrix const & | params = base::DataMatrix(0, 0) |
||
) |
See the other version of evaluate() and the setParameters() overloads.
std::shared_ptr< AbstractMultiStorage< FloatArrayVector > > sgpp::combigrid::CombigridMultiOperation::getDifferences | ( | ) |
Returns a storage of the differences (Deltas) that have been computed by the CombigridEvaluator.
std::shared_ptr< LevelManager > sgpp::combigrid::CombigridMultiOperation::getLevelManager | ( | ) |
Via the LevelManager, more options are available than are provided directly by this class.
std::vector< std::shared_ptr< AbstractPointHierarchy > > sgpp::combigrid::CombigridMultiOperation::getPointHierarchies | ( | ) |
base::DataVector sgpp::combigrid::CombigridMultiOperation::getResult | ( | ) |
Returns the current computed value.
References python.statsfileInfo::i.
Referenced by python.uq.uq_setting.UQSetting.UQSetting::getTimeDependentResults().
std::shared_ptr< AbstractCombigridStorage > sgpp::combigrid::CombigridMultiOperation::getStorage | ( | ) |
size_t sgpp::combigrid::CombigridMultiOperation::getUpperPointBound | ( | ) | const |
size_t sgpp::combigrid::CombigridMultiOperation::numDims | ( | ) |
size_t sgpp::combigrid::CombigridMultiOperation::numGridPoints | ( | ) |
size_t sgpp::combigrid::CombigridMultiOperation::numStoredFunctionValues | ( | ) |
void sgpp::combigrid::CombigridMultiOperation::setLevelManager | ( | std::shared_ptr< LevelManager > | levelManager | ) |
Can be used to set the level manager, e.g.
if one of the static constructor functions has been used.
References sgpp::combigrid::CombigridMultiOperationImpl::levelManager.
void sgpp::combigrid::CombigridMultiOperation::setParameters | ( | std::vector< base::DataVector > const & | params | ) |
Sets the parameters for upcoming computations and clears the data structures (removes old computed data).
This is only relevant for methods with parameters (e.g. interpolation, but not quadrature).
params | The parameters at which the function should be evaluated. |
References python.statsfileInfo::i, and python.utils.statsfile2gnuplot::j.
void sgpp::combigrid::CombigridMultiOperation::setParameters | ( | base::DataMatrix const & | params | ) |
Sets the parameters for upcoming computations and clears the data structures (removes old computed data).
This is only relevant for methods with parameters (e.g. interpolation, but not quadrature). The evaluation points are assumed to be the columns of the matrix.
References sgpp::base::DataMatrix::getNcols(), sgpp::base::DataMatrix::getNrows(), python.statsfileInfo::i, and python.utils.statsfile2gnuplot::j.