31 #ifndef __TSG_CACHE_LAGRANGE_HPP
32 #define __TSG_CACHE_LAGRANGE_HPP
76 cache(std::vector<std::vector<T>>(num_dimensions, std::vector<T>())),
offsets(rule.getPointsCount()){
77 for(
int dim=0; dim<num_dimensions; dim++){
79 for(
int level=0; level <= max_levels[dim]; level++)
88 const double *nodes = rule.
getNodes(level);
94 for(
int j=0; j<num_points-1; j++){
99 cc[num_points-1] *= c * coeff[num_points-1];
100 for(
int j=num_points-2; j>=0; j--){
101 c *= (x - nodes[j+1]);
102 cc[j] *= c * coeff[j];
133 template <
typename T>
149 cache(std::vector<std::vector<T>>(num_dimensions, std::vector<T>())),
offsets(rule.getPointsCount()) {
150 for(
int dim=0; dim<num_dimensions; dim++){
152 for(
int level=0; level <= max_levels[dim]; level++)
170 const double *nodes = rule.
getNodes(level);
174 std::vector<T> aux_f(num_points), aux_g(num_points);
177 aux_g[num_points-1] = 1.0;
178 for(
int i=1; i<num_points; i++) {
179 aux_f[i] = aux_f[i-1] * (x - nodes[i-1]);
180 aux_g[num_points-1-i] = aux_g[num_points-i] * (x - nodes[num_points-i]);
181 cc[i] = aux_f[i-1] + (x - nodes[i-1]) * cc[i-1];
183 cc[num_points-1] *= coeff[num_points-1];
185 for(
int i=num_points-2; i>=0; i--) {
186 diff_gj = aux_g[i+1] + diff_gj * (x - nodes[i+1]);
187 cc[i] = coeff[i] * (cc[i] * aux_g[i] + aux_f[i] * diff_gj);
Cache that holds the derivatives of 1D Lagrange polynomials. Uses the same interface as CacheLagrange...
Definition: tsgCacheLagrange.hpp:134
std::vector< int > offsets
See CacheLagrange::offsets.
Definition: tsgCacheLagrange.hpp:200
static void cacheDerivativeLevel(int level, double x, const OneDimensionalWrapper &rule, T *cc)
Computes the derivatives of all Lagrange polynomials for the given level at the given x.
Definition: tsgCacheLagrange.hpp:160
~CacheLagrangeDerivative()=default
Destructor, clear all used data.
T getLagrangeDerivative(int dimension, int level, int local) const
Return the Lagrange derivative cache for given dimension, level and offset local to the level.
Definition: tsgCacheLagrange.hpp:192
CacheLagrangeDerivative(int num_dimensions, const std::vector< int > &max_levels, const OneDimensionalWrapper &rule, const double x[])
Constructor that takes into account a single canonical point x.
Definition: tsgCacheLagrange.hpp:148
std::vector< std::vector< T > > cache
See CacheLagrange::cache.
Definition: tsgCacheLagrange.hpp:198
Cache that holds the values of 1D Lagrange polynomials.
Definition: tsgCacheLagrange.hpp:61
std::vector< std::vector< T > > cache
Stores the values of the Lagrange polynomias for each dimension and each level.
Definition: tsgCacheLagrange.hpp:120
std::vector< int > offsets
Stores the offsets for each level within each dimension.
Definition: tsgCacheLagrange.hpp:122
static void cacheLevel(int level, double x, const OneDimensionalWrapper &rule, T *cc)
Computes the values of all Lagrange polynomials for the given level at the given x.
Definition: tsgCacheLagrange.hpp:87
T getLagrange(int dimension, int level, int local) const
Return the Lagrange cache for given dimension, level and offset local to the level.
Definition: tsgCacheLagrange.hpp:107
~CacheLagrange()=default
Destructor, clear all used data.
CacheLagrange(int num_dimensions, const std::vector< int > &max_levels, const OneDimensionalWrapper &rule, const double x[])
Constructor that takes into account a single canonical point x.
Definition: tsgCacheLagrange.hpp:75
A class to cache one dimensional rules, nodes, weight, meta-data, etc.
Definition: tsgOneDimensionalWrapper.hpp:53
const double * getCoefficients(int level) const
Get an array of the Lagrange coefficients associated with the level, the order matches getNodes().
Definition: tsgOneDimensionalWrapper.hpp:305
const double * getNodes(int level) const
Get an array of all nodes associated with the level (the array is ordered by local index).
Definition: tsgOneDimensionalWrapper.hpp:303
TypeOneDRule getRule() const
Get the loaded rule.
Definition: tsgOneDimensionalWrapper.hpp:311
int getNumPoints(int level) const
Get the number of points for the level, can only query loaded levels.
Definition: tsgOneDimensionalWrapper.hpp:279
@ rule_clenshawcurtis0
Same as rule_clenshawcurtis but with modified basis that assumes the model is zero at the boundary.
Definition: tsgEnumerates.hpp:291
Encapsulates the Tasmanian Sparse Grid module.
Definition: TasmanianSparseGrid.hpp:68
Definition for the base grid class.