Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2 (development)
TasGrid::OneDimensionalNodes Namespace Reference

Contains algorithms for generating Gauss,Chebyshev, and Fourier nodes and weights. More...

Functions

void getGaussLegendre (int m, std::vector< double > &w, std::vector< double > &x)
 Generate Gauss-Legendre weights w and points x for (input) number of points m.
 
void getChebyshev (int m, std::vector< double > &w, std::vector< double > &x)
 Generate Chebyshev weights w and points x for (input) number of points m.
 
void getGaussChebyshev1 (int m, std::vector< double > &w, std::vector< double > &x)
 Generate Gauss-Chebyshev type 1 weights w and points x for (input) number of points m.
 
void getGaussChebyshev2 (int m, std::vector< double > &w, std::vector< double > &x)
 Generate Gauss-Chebyshev type 2 weights w and points x for (input) number of points m.
 
void getGaussJacobi (int m, std::vector< double > &w, std::vector< double > &x, double alpha, double beta)
 Generate Gauss-Jacobi weights w and points x for (input) number of points m, using parameters alpha and beta.
 
void getGaussHermite (int m, std::vector< double > &w, std::vector< double > &x, double alpha)
 Generate Gauss-Hermite weights w and points x for (input) number of points m, using parameters alpha.
 
void getGaussLaguerre (int m, std::vector< double > &w, std::vector< double > &x, double alpha)
 Generate Gauss-Laguerre weights w and points x for (input) number of points m, using parameters alpha.
 
std::vector< double > getClenshawCurtisNodes (int level)
 Generate Clenshaw-Curtis nodes for the level.
 
double getClenshawCurtisWeight (int level, int point)
 Return the Clenshaw-Curtis weight for the level and node indexed by point.
 
std::vector< double > getClenshawCurtisNodesZero (int level)
 Generate Clenshaw-Curtis zero-boundary condition nodes for the level.
 
double getClenshawCurtisWeightZero (int level, int point)
 Return the Clenshaw-Curtis zero-boundary condition weight for the level and node indexed by point.
 
std::vector< double > getFejer2Nodes (int level)
 Generate Fejer type 2 nodes for the level.
 
double getFejer2Weight (int level, int point)
 Return the Fejer type 2 weight for the level and node indexed by point.
 
std::vector< double > getRLeja (int n)
 Generate the first n R-Leja nodes, starting with 1, -1, 0, ...
 
std::vector< double > getRLejaCentered (int n)
 Generate the first n R-Leja nodes, starting with 0, 1, -1, ...
 
std::vector< double > getRLejaShifted (int n)
 Generate the first n R-Leja nodes, starting with -0.5, 0.5, ...
 
std::vector< double > getFourierNodes (int level)
 Generate the Fourier nodes for the given level, uniformly distributed points with right-most point omitted due to periodicity.
 

Detailed Description

Contains algorithms for generating Gauss,Chebyshev, and Fourier nodes and weights.