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

Files

file  tsgExoticQuadrature.hpp
 Header to include the exotic quadrature templates.
 

Functions

double TasGrid::poly_eval (const std::vector< double > &roots, double x)
 Evaluates a polynomial with roots given by roots at the point x. More...
 
double TasGrid::lagrange_eval (size_t idx, const std::vector< double > &roots, double x)
 Evaluates a Lagrange basis polynomial at a point x. More...
 
template<bool is_symmetric>
void TasGrid::getGaussNodesAndWeights (const int n, const std::vector< double > &ref_points, const std::vector< double > &ref_weights, std::vector< std::vector< double >> &points_cache, std::vector< std::vector< double >> &weights_cache)
 Generates n levels of points and weights using orthogonal polynomials with respect to ref_points and ref_weights. More...
 
TasGrid::CustomTabulated TasGrid::getShiftedExoticQuadrature (const int n, const double shift, const std::vector< double > &shifted_weights, const std::vector< double > &ref_points, const char *description, const bool is_symmetric=false)
 Generate the exotic quadrature points and weights and load them into a TasGrid::CustomTabulated object. More...
 
void TasGrid::shiftReferenceWeights (std::vector< double > const &vals, double shift, std::vector< double > &ref_weights)
 Multiplies the reference quadrature weights by the shifted values of the weight function. More...
 
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

◆ poly_eval()

double TasGrid::poly_eval ( const std::vector< double > &  roots,
double  x 
)
inline

Evaluates a polynomial with roots given by roots at the point x.

The polynomial is defined by a simple product of (x - roots[i]) and will not be normalized.

◆ lagrange_eval()

double TasGrid::lagrange_eval ( size_t  idx,
const std::vector< double > &  roots,
double  x 
)
inline

Evaluates a Lagrange basis polynomial at a point x.

The Lagrange basis at idx is 1 when x = roots[idx] and 0 for all other roots.

Parameters
idxis the index of the Lagrange polynomial
rootsis the point set the defines the Lagrange polynomials
xis the point of interest

◆ getGaussNodesAndWeights()

template<bool is_symmetric>
void TasGrid::getGaussNodesAndWeights ( const int  n,
const std::vector< double > &  ref_points,
const std::vector< double > &  ref_weights,
std::vector< std::vector< double >> &  points_cache,
std::vector< std::vector< double >> &  weights_cache 
)
inline

Generates n levels of points and weights using orthogonal polynomials with respect to ref_points and ref_weights.

Given the reference points and weights, finds the roots of the first n + 1 orthogonal polynomials. The first polynomial is always constant, every following polynomial is chosen so that it is orthogonal with respect to all previous polynomials and the orthogonality (i.e., inner product integral) is computed with the reference quadrature points and weights.

Template Parameters
is_symmetricindicates whether the weight function is symmetric or not, the algorithm makes slight modifications that improve stability in the symmetric case
Parameters
nis the number of quadrature rules to construct
ref_pointsare the points of the reference quadrature
ref_weightsare the weights of the reference quadrature multiplied by the shifted weight function
points_cacheon exit will be loaded so that points_cache[i] will be roots of the i+1-st orthogonal polynomial
weights_cacheon exit will be loaded with the quadrature weights corresponding to the points_cache

◆ getShiftedExoticQuadrature()

TasGrid::CustomTabulated TasGrid::getShiftedExoticQuadrature ( const int  n,
const double  shift,
const std::vector< double > &  shifted_weights,
const std::vector< double > &  ref_points,
const char *  description,
const bool  is_symmetric = false 
)
inline

Generate the exotic quadrature points and weights and load them into a TasGrid::CustomTabulated object.

The method calls getGaussNodesAndWeights() and adds the Gauss-Legendre correction.

Parameters
nis the number of levels to compute
shiftis the correction term applied to the weight function to enforce positivity
shifted_weightsare the reference quadrature points multiplied by the values of the shifted weight function at the reference points
ref_pointsare the reference quadrature points
descriptionis the human readable string to identify the quadrature rule
is_symmetricindicates whether we can use algorithm modifications that can improve stability

◆ shiftReferenceWeights()

void TasGrid::shiftReferenceWeights ( std::vector< double > const &  vals,
double  shift,
std::vector< double > &  ref_weights 
)
inline

Multiplies the reference quadrature weights by the shifted values of the weight function.

Parameters
valsare the values of the exotic weight function at the reference quadrature points
shiftis a real-valued scalar that ensures all vals[i] + shift are positive
ref_weightsare the reference quadrature weights, will be overwritten with the result ref_weights[i] * (vals[i] + shift)
Exceptions
std::invalid_argumentif vals[i] + shift is negative for some index i

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