|
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...
|
|
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_symmetric | indicates whether the weight function is symmetric or not, the algorithm makes slight modifications that improve stability in the symmetric case |
- Parameters
-
n | is the number of quadrature rules to construct |
ref_points | are the points of the reference quadrature |
ref_weights | are the weights of the reference quadrature multiplied by the shifted weight function |
points_cache | on exit will be loaded so that points_cache[i] will be roots of the i+1-st orthogonal polynomial |
weights_cache | on exit will be loaded with the quadrature weights corresponding to the points_cache |
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_levels | is the number of levels that should be computed |
shift | is a real-valued scalar that needs to be added to the weight function to make it positive |
grid | has 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) |
description | is a human readable string that can be used to identify the quadrature rule (could be empty) |
is_symmetric | indicates 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_argument | if 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
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
,
, constant weight, and the Gauss-Chebyshev weights. On the other hand,
is not symmetric in the context of this algorithm due to the shift that must be applied.