SG++-Doxygen-Documentation
sgpp::combigrid::CombigridMultiOperation Class Reference

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 &params)
 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 &params=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< LevelManagergetLevelManager ()
 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< AbstractCombigridStoragegetStorage ()
 
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 &params)
 Sets the parameters for upcoming computations and clears the data structures (removes old computed data). More...
 
void setParameters (base::DataMatrix const &params)
 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::CombigridMultiOperationcreateBsplineLinearCoefficientOperation (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::CombigridMultiOperationcreateBsplineLinearRefinementOperation (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::CombigridMultiOperationcreateBsplineVarianceRefinementOperation (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< CombigridMultiOperationcreateExpChebyshevPolynomialInterpolation (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< CombigridMultiOperationcreateExpClenshawCurtisPolynomialInterpolation (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< CombigridMultiOperationcreateExpClenshawCurtisQuadrature (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< CombigridMultiOperationcreateExpL2LejaPolynomialInterpolation (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< CombigridMultiOperationcreateExpLejaPolynomialInterpolation (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< CombigridMultiOperationcreateExpUniformBoundaryBsplineInterpolation (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< CombigridMultiOperationcreateExpUniformBoundaryBsplineQuadrature (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< CombigridMultiOperationcreateExpUniformBoundaryBsplineSquareQuadrature (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< CombigridMultiOperationcreateExpUniformBoundaryLinearInterpolation (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< CombigridMultiOperationcreateExpUniformLinearInterpolation (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< CombigridMultiOperationcreateLinearClenshawCurtisPolynomialInterpolation (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< CombigridMultiOperationcreateLinearL2LejaPolynomialInterpolation (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< CombigridMultiOperationcreateLinearL2LejaQuadrature (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< CombigridMultiOperationcreateLinearLejaPolynomialInterpolation (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< CombigridMultiOperationcreateLinearLejaQuadrature (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...
 

Detailed Description

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:

  • The point hierarchies and evaluators etc. are created by the user and passed to the constructor
  • One of the static methods is used. They provide some sensible isotropic configurations.

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:

Constructor & Destructor Documentation

◆ CombigridMultiOperation() [1/3]

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.

◆ CombigridMultiOperation() [2/3]

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.

◆ CombigridMultiOperation() [3/3]

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.

Member Function Documentation

◆ createBsplineLinearCoefficientOperation()

std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createBsplineLinearCoefficientOperation ( size_t  degree,
size_t  numDimensions,
std::shared_ptr< sgpp::combigrid::AbstractCombigridStorage coefficientStorage 
)
static

creates a B spline interpolation operation from a storage of interpolation coefficients

Parameters
degreedegree of the b-splines
numDimensionsnumber of dimensions
coefficientStoragestorage of the b-spline coefficients
Returns
combigrid multi operation

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().

◆ createBsplineLinearRefinementOperation()

std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createBsplineLinearRefinementOperation ( size_t  degree,
size_t  numDimensions,
sgpp::combigrid::MultiFunction  func,
std::shared_ptr< sgpp::combigrid::LevelManager levelManager 
)
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

Parameters
degreeB spline degree
numDimensionsnumber of dimensions
functhe objective function
levelManagerlevel 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.

◆ createBsplineVarianceRefinementOperation()

std::shared_ptr< sgpp::combigrid::CombigridMultiOperation > 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 
)
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

Parameters
degreeB spline degree
numDimensionsnumber of dimensions
functhe objective function
levelManagerlevel manager
weightFunctionsweight functions
boundsbounding box
Returns
a combigrid operation calculating the variance on each full grid

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.

◆ createExpChebyshevPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpChebyshevPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a Chebyshev grid with an exponential growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expChebyshev(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createExpClenshawCurtisPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpClenshawCurtisPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with an exponential growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expClenshawCurtis(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createExpClenshawCurtisQuadrature()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpClenshawCurtisQuadrature ( size_t  numDimensions,
MultiFunction  func 
)
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.

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be integrated.

References sgpp::combigrid::CombiHierarchies::expClenshawCurtis(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().

◆ createExpL2LejaPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpL2LejaPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with an exponential growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expL2Leja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createExpLejaPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpLejaPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with an exponential growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expLeja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createExpUniformBoundaryBsplineInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpUniformBoundaryBsplineInterpolation ( size_t  numDimensions,
MultiFunction  func,
size_t  degree = 3 
)
static

Returns a CombigridMultiOperation doing B-spline interpolation on a uniform grid with boundary using an exponential growth strategy (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.
degreedegree 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.

◆ createExpUniformBoundaryBsplineQuadrature()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpUniformBoundaryBsplineQuadrature ( size_t  numDimensions,
MultiFunction  func,
size_t  degree 
)
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.

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be integrated.
degreeB-spline degree

References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombiEvaluators::multiBSplineQuadrature(), and sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies.

◆ createExpUniformBoundaryBsplineSquareQuadrature()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpUniformBoundaryBsplineSquareQuadrature ( size_t  numDimensions,
MultiFunction  func,
size_t  degree 
)
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.

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be integrated.
degreeB-spline degree

References BSplineCoefficientGridFunction(), sgpp::combigrid::CombiEvaluators::BSplineScalarProduct(), sgpp::combigrid::CombiHierarchies::expUniformBoundary(), sgpp::combigrid::CombigridMultiOperationImpl::levelManager, sgpp::combigrid::CombigridMultiOperationImpl::pointHierarchies, and sgpp::combigrid::QUADRATIC.

◆ createExpUniformBoundaryLinearInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpUniformBoundaryLinearInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid with boundary using an exponential growth strategy (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expUniformBoundary(), and sgpp::combigrid::CombiEvaluators::multiLinearInterpolation().

◆ createExpUniformLinearInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createExpUniformLinearInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing (multi-)linear interpolation on a uniform grid without boundary using an exponential growth strategy (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiHierarchies::expUniform(), and sgpp::combigrid::CombiEvaluators::multiLinearInterpolation().

◆ createLinearClenshawCurtisPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createLinearClenshawCurtisPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a Clenshaw-Curtis grid with a linear growth (not nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.

References sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createLinearL2LejaPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createLinearL2LejaPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func,
size_t  growthFactor = 2 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a L2Leja grid with linear growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.
growthFactorParameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used.

References sgpp::combigrid::CombiHierarchies::linearL2Leja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createLinearL2LejaQuadrature()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createLinearL2LejaQuadrature ( size_t  numDimensions,
MultiFunction  func,
size_t  growthFactor = 2 
)
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.

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be integrated.
growthFactorParameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used.

References sgpp::combigrid::CombiHierarchies::linearL2Leja(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().

◆ createLinearLejaPolynomialInterpolation()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createLinearLejaPolynomialInterpolation ( size_t  numDimensions,
MultiFunction  func,
size_t  growthFactor = 2 
)
static

Returns a CombigridMultiOperation doing polynomial interpolation on a Leja grid with linear growth (nested points).

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be interpolated.
growthFactorParameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used.

References sgpp::combigrid::CombiHierarchies::linearLeja(), and sgpp::combigrid::CombiEvaluators::multiPolynomialInterpolation().

◆ createLinearLejaQuadrature()

std::shared_ptr< CombigridMultiOperation > sgpp::combigrid::CombigridMultiOperation::createLinearLejaQuadrature ( size_t  numDimensions,
MultiFunction  func,
size_t  growthFactor = 2 
)
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.

Parameters
numDimensionsDimensionality of the problem.
funcFunction to be integrated.
growthFactorParameter for the linear growth strategy. For level l, 1 + growthFactor * l points are used.

References sgpp::combigrid::CombiHierarchies::linearLeja(), and sgpp::combigrid::CombiEvaluators::multiQuadrature().

◆ evaluate() [1/2]

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.

◆ evaluate() [2/2]

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.

◆ getDifferences()

std::shared_ptr< AbstractMultiStorage< FloatArrayVector > > sgpp::combigrid::CombigridMultiOperation::getDifferences ( )

Returns a storage of the differences (Deltas) that have been computed by the CombigridEvaluator.

◆ getLevelManager()

std::shared_ptr< LevelManager > sgpp::combigrid::CombigridMultiOperation::getLevelManager ( )

Via the LevelManager, more options are available than are provided directly by this class.

◆ getPointHierarchies()

std::vector< std::shared_ptr< AbstractPointHierarchy > > sgpp::combigrid::CombigridMultiOperation::getPointHierarchies ( )
Returns
the point hierarchies containing the grid points in each direction

◆ getResult()

base::DataVector sgpp::combigrid::CombigridMultiOperation::getResult ( )

Returns the current computed value.

References python.statsfileInfo::i.

Referenced by python.uq.uq_setting.UQSetting.UQSetting::getTimeDependentResults().

◆ getStorage()

std::shared_ptr< AbstractCombigridStorage > sgpp::combigrid::CombigridMultiOperation::getStorage ( )
Returns
the storage containing the computed coefficients. For the basic operations these are the function values at evaluation points.

◆ getUpperPointBound()

size_t sgpp::combigrid::CombigridMultiOperation::getUpperPointBound ( ) const
Returns
An upper bound for the number of points (function evaluations) used for the current computation. This bound is exact if nesting is used or if otherwise each grid point only occurs in exactly one level.

◆ numDims()

size_t sgpp::combigrid::CombigridMultiOperation::numDims ( )
Returns
the number of dimensions

◆ numGridPoints()

size_t sgpp::combigrid::CombigridMultiOperation::numGridPoints ( )
Returns
the total number of different (multi-dimensional) grid points that have been used for the current evaluation. This number is reset when clear() is called on the CombigridEvaluator via evaluate() or setParameters(). This method is currently not optimized and can be slow!

◆ numStoredFunctionValues()

size_t sgpp::combigrid::CombigridMultiOperation::numStoredFunctionValues ( )
Returns
the number of function values that have been computed via this CombigridMultiOperation during its lifetime. For a nested grid, this number matches numGridPoints() if only one computation is performed, i.e. no previous data has been cleared via evaluate() or setParameters(). Its computation is not optimized, but currently faster than numGridPoints().

◆ setLevelManager()

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.

◆ setParameters() [1/2]

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

Parameters
paramsThe parameters at which the function should be evaluated.

References python.statsfileInfo::i, and python.utils.statsfile2gnuplot::j.

◆ setParameters() [2/2]

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.


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