Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.1
Exotic Quadrature
Collaboration diagram for Exotic Quadrature:

Functions

TasGrid::CustomTabulated TasGrid::getExoticQuadrature (const int num_levels, const double shift, const TasGrid::TasmanianSparseGrid &grid, const char *description, const bool is_symmetric=false)
 Constructs an exotic quadrature from a sparse grid object. More...
 
TasGrid::CustomTabulated TasGrid::getExoticQuadrature (const int num_levels, const double shift, std::function< double(double)> weight_fn, const int num_ref_points, const char *description, const bool is_symmetric=false)
 Similar to getExoticQuadrature() but the weight function is defined by a lambda expression. More...
 

Detailed Description

Exotic Quadrature
The Exotic Quadrature module of Tasmanian offers a collection of functions for generating one dimensional quadrature rules with weight functions that are not strictly non-negative. The quadrature rules are loaded inside TasGrid::CustomTabulated objects and be used in calls to TasmanianSparseGrid::makeGlobalGrid() to extend those to a multidimensional context.
The integral of interest has the form $ \int_{-1}^{1} f(x) \rho(x) dx $ and it differs from the classical Gaussian rules by allowing the weight function $ \rho(x) $ to be negative in portions of the domain. We still require that the weight is integrable and bounded from below, i.e., $ \int_{-1}^{1} | \rho(x) | dx < \infty $ and $ \rho(x) \geq c > -\infty $.
The quadrature consists of two components, one is a quadrature that uses classical Gaussian construction using the root of polynomials that are orthogonal with respect to a shifted positive weight $ \rho(x) + c \geq 0 $ and a correction term using Gauss-Legendre rule scaled by -c. While the asymptotic convergence rate is inferior to that of Gauss-Legendre, for small number of points the accuracy of the combined quadrature rule can be significantly greater. This advantage is further magnified in a sparse grid context where reaching the asymptotic regime for Gauss-Legendre may require a prohibitive number of points due to the curse of dimensionality.
The weight can be defined either with a C++ lambda expression or a one-dimensional TasmanianSparseGrid object. The lambda expression will be evaluated at a large number of Gauss-Legendre points that will be used for a reference quadrature to build the orthogonal polynomial basis. In the case of a sparse grid object, only the loaded points and quadrature weights will be used, and no interpolation will be performed inbetween (i.e., we circumvent potential interpolation error).

Function Documentation

◆ getExoticQuadrature() [1/2]

TasGrid::CustomTabulated TasGrid::getExoticQuadrature ( const int  num_levels,
const double  shift,
const TasGrid::TasmanianSparseGrid grid,
const char *  description,
const bool  is_symmetric = false 
)
inline

Constructs an exotic quadrature from a sparse grid object.

Constructs a set of one dimensional quadrature rules with respect to a 1D weight function approximated by grid.

Parameters
num_levelsis the number of levels that should be computed
shiftis a real-valued scalar that needs to be added to the weight function to make it positive
gridhas dimension 1 and at least one loaded point, should use as many points as possible for better accuracy (up to some point where numerical stability is lost)
descriptionis a human readable string that can be used to identify the quadrature rule (could be empty)
is_symmetricindicates whether to assume that the weight is symmetric, which allows for some stability improvements; symmetric means "even" on [-1, 1], which leads to all odd power polynomials integrating to zero; the same holds for shifted domains so long as the weight is "even" with respect to the mid-point of the domain.
Returns
TasGrid::CustomTabulated object holding the points and weights for the different levels of the quadrature
Exceptions
std::invalid_argumentif grid doesn't have a single input and output or no loaded points, or if the value of the weight function plus the shift is not non-negative

The grid object is used in two ways: first, it defines a reference quadrature that will be used to compute the inner-product integrals that define orthogonal polynomials; second, it defines the values of the weight function $ \rho(x) $ at the quadrature points. Note that the values of the weight function will not be interpolated between the reference points.

Examples of symmetric weight functions are $ \sin(x) $, $ \sin(x) / x $, constant weight, and the Gauss-Chebyshev weights. On the other hand, $ \cos(x) $ is not symmetric in the context of this algorithm due to the shift that must be applied.

◆ getExoticQuadrature() [2/2]

TasGrid::CustomTabulated TasGrid::getExoticQuadrature ( const int  num_levels,
const double  shift,
std::function< double(double)>  weight_fn,
const int  num_ref_points,
const char *  description,
const bool  is_symmetric = false 
)
inline

Similar to getExoticQuadrature() but the weight function is defined by a lambda expression.

Constructs a CustomTabulated object holding the points and weights for the different levels of the quadrature for the given parameter weight_fn. The num_ref_points parameter sets the number of reference points from a Gauss-Legendre quadrature that will be used to compute the inner-product integrals that define the orthogonal polynomials.