Doxygen
1.9.1
|
Algorithms and meta-data for the Gauuss, Chebyshev and custom rules. More...
Go to the source code of this file.
Classes | |
class | TasGrid::CustomTabulated |
Class providing manipulation of custom tabulated rules, file I/O and structured access to the points, weights, and meta-data. More... | |
Namespaces | |
TasGrid | |
Encapsulates the Tasmanian Sparse Grid module. | |
TasGrid::OneDimensionalMeta | |
Metadata for one dimensional rules, number of points, exactness, I/O, etc. | |
TasGrid::OneDimensionalNodes | |
Contains algorithms for generating Gauss,Chebyshev, and Fourier nodes and weights. | |
Functions | |
CustomTabulated | TasGrid::getSubrules (CustomTabulated &ct, int start_index, int stride, std::string description) |
Generates a subset of rules as a CustomTabulated object. The subset has a short description string given by description and takes the levels of the input ct starting from level start_index with distance stride between consecutive levels. More... | |
int | TasGrid::OneDimensionalMeta::getNumPoints (int level, TypeOneDRule rule) |
Return the number of points for the rule at the level, includes all global rules. | |
int | TasGrid::OneDimensionalMeta::getIExact (int level, TypeOneDRule rule) |
Return the exactness of the interpolation rule at the level, includes all global rules. | |
int | TasGrid::OneDimensionalMeta::getQExact (int level, TypeOneDRule rule) |
Return the exactness of the integration/quadrature rule at the level, includes all global rules. | |
bool | TasGrid::OneDimensionalMeta::isNonNested (TypeOneDRule rule) |
Return True if the rule does not have nested nodes, e.g., gauss-legendre. | |
bool | TasGrid::OneDimensionalMeta::isSequence (TypeOneDRule rule) |
Return True if the rule is nested and has single-node growth, fit for GridSequence, e.g., leja. | |
bool | TasGrid::OneDimensionalMeta::isGlobal (TypeOneDRule rule) |
Return True if the rule has basis with Lagrange polynomials with global support, fit for GridGlobal, e.g., clenshaw-curtis. | |
bool | TasGrid::OneDimensionalMeta::isSingleNodeGrowth (TypeOneDRule rule) |
Return True if the rule grows by one point per level. | |
bool | TasGrid::OneDimensionalMeta::isLocalPolynomial (TypeOneDRule rule) |
Return True if the rule has polynomial basis with local support, fit for GridLocalPolynomial, e.g., localp. | |
bool | TasGrid::OneDimensionalMeta::isWavelet (TypeOneDRule rule) |
Return True if the rule has wavelet basis with local support. | |
bool | TasGrid::OneDimensionalMeta::isFourier (TypeOneDRule rule) |
Return True if the rule has trigonometric basis. | |
const char * | TasGrid::OneDimensionalMeta::getHumanString (TypeOneDRule rule) |
Map the enumerate to a human readable string, used in printStats(). | |
TypeDepth | TasGrid::OneDimensionalMeta::getControurType (TypeDepth type) |
Identifies the general contour, linear type_level, log-corrected type_curved, or hyperbolic type_hyperbolic. More... | |
bool | TasGrid::OneDimensionalMeta::isExactLevel (TypeDepth type) |
Returns true if the type indicates exactness with respect to raw levels. | |
bool | TasGrid::OneDimensionalMeta::isExactInterpolation (TypeDepth type) |
Returns true if the type indicates exactness with respect to interpolation. | |
bool | TasGrid::OneDimensionalMeta::isExactQuadrature (TypeDepth type) |
Returns true if the type indicates exactness with respect to integration. | |
TypeDepth | TasGrid::OneDimensionalMeta::getSelectionType (TypeDepth type) |
Identifies the selection type, level type_level, interpolation type_iptotal, or quadrature type_qptotal. More... | |
bool | TasGrid::OneDimensionalMeta::isTypeCurved (TypeDepth type) |
Return True if the multi-index selection type has log-correction term (need to use floating point indexing). | |
void | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 | TasGrid::OneDimensionalNodes::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 > | TasGrid::OneDimensionalNodes::getClenshawCurtisNodes (int level) |
Generate Clenshaw-Curtis nodes for the level. | |
double | TasGrid::OneDimensionalNodes::getClenshawCurtisWeight (int level, int point) |
Return the Clenshaw-Curtis weight for the level and node indexed by point. | |
std::vector< double > | TasGrid::OneDimensionalNodes::getClenshawCurtisNodesZero (int level) |
Generate Clenshaw-Curtis zero-boundary condition nodes for the level. | |
double | TasGrid::OneDimensionalNodes::getClenshawCurtisWeightZero (int level, int point) |
Return the Clenshaw-Curtis zero-boundary condition weight for the level and node indexed by point. | |
std::vector< double > | TasGrid::OneDimensionalNodes::getFejer2Nodes (int level) |
Generate Fejer type 2 nodes for the level. | |
double | TasGrid::OneDimensionalNodes::getFejer2Weight (int level, int point) |
Return the Fejer type 2 weight for the level and node indexed by point. | |
std::vector< double > | TasGrid::OneDimensionalNodes::getRLeja (int n) |
Generate the first n R-Leja nodes, starting with 1, -1, 0, ... | |
std::vector< double > | TasGrid::OneDimensionalNodes::getRLejaCentered (int n) |
Generate the first n R-Leja nodes, starting with 0, 1, -1, ... | |
std::vector< double > | TasGrid::OneDimensionalNodes::getRLejaShifted (int n) |
Generate the first n R-Leja nodes, starting with -0.5, 0.5, ... | |
std::vector< double > | TasGrid::OneDimensionalNodes::getFourierNodes (int level) |
Generate the Fourier nodes for the given level, uniformly distributed points with right-most point omitted due to periodicity. | |
Algorithms and meta-data for the Gauuss, Chebyshev and custom rules.
Contains core information about one dimensional rules, custom tabulated, Chebyshev and Gaussian rules. Also, generic I/O data.