Stokhos
Development
|
Abstract base class for multivariate orthogonal polynomials. More...
#include <Stokhos_OrthogPolyBasis.hpp>
Public Member Functions | |
OrthogPolyBasis () | |
Constructor. | |
virtual | ~OrthogPolyBasis () |
Destructor. | |
virtual ordinal_type | order () const =0 |
Return order of basis. | |
virtual ordinal_type | dimension () const =0 |
Return dimension of basis. | |
virtual ordinal_type | size () const =0 |
Return total size of basis. | |
virtual const Teuchos::Array < value_type > & | norm_squared () const =0 |
Return array storing norm-squared of each basis polynomial. More... | |
virtual const value_type & | norm_squared (ordinal_type i) const =0 |
Return norm squared of basis polynomial i . | |
virtual Teuchos::RCP < Stokhos::Sparse3Tensor < ordinal_type, value_type > > | computeTripleProductTensor () const =0 |
Compute triple product tensor. More... | |
virtual Teuchos::RCP < Stokhos::Sparse3Tensor < ordinal_type, value_type > > | computeLinearTripleProductTensor () const =0 |
Compute linear triple product tensor where k = 0,1. | |
virtual value_type | evaluateZero (ordinal_type i) const =0 |
Evaluate basis polynomial i at zero. | |
virtual void | evaluateBases (const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const =0 |
Evaluate basis polynomials at given point point . More... | |
virtual void | print (std::ostream &os) const =0 |
Print basis to stream os . | |
virtual const std::string & | getName () const =0 |
Return string name of basis. | |
Abstract base class for multivariate orthogonal polynomials.
This class provides an abstract interface for multivariate orthogonal polynomials. Orthogonality is defined by the inner product
where is the density function of the measure associated with the orthogonal polynomials and is the dimension of the domain.
Like most classes in Stokhos, the class is templated on the ordinal and value types. Typically ordinal_type
= int
and value_type
= double
.
|
pure virtual |
Compute triple product tensor.
The entry of the tensor is given by where represents basis polynomial and where is size()-1.
Implemented in Stokhos::CompletePolynomialBasis< ordinal_type, value_type >, Stokhos::ProductLanczosPCEBasis< ordinal_type, value_type >, Stokhos::GSReducedPCEBasisBase< ordinal_type, value_type >, Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >, Stokhos::ProductLanczosGramSchmidtPCEBasis< ordinal_type, value_type >, Stokhos::TensorProductBasis< ordinal_type, value_type, coeff_compare_type >, Stokhos::TotalOrderBasis< ordinal_type, value_type, coeff_compare_type >, and Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >.
|
pure virtual |
Evaluate basis polynomials at given point point
.
Size of returned array is given by size(), and coefficients are ordered from order 0 up to size size()-1.
Implemented in Stokhos::CompletePolynomialBasis< ordinal_type, value_type >, Stokhos::ProductLanczosPCEBasis< ordinal_type, value_type >, Stokhos::GSReducedPCEBasisBase< ordinal_type, value_type >, Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >, Stokhos::ProductLanczosGramSchmidtPCEBasis< ordinal_type, value_type >, Stokhos::TensorProductBasis< ordinal_type, value_type, coeff_compare_type >, Stokhos::TotalOrderBasis< ordinal_type, value_type, coeff_compare_type >, and Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >.
|
pure virtual |
Return array storing norm-squared of each basis polynomial.
Entry of returned array is given by for where is size()-1.
Implemented in Stokhos::CompletePolynomialBasis< ordinal_type, value_type >, Stokhos::ProductLanczosPCEBasis< ordinal_type, value_type >, Stokhos::GSReducedPCEBasisBase< ordinal_type, value_type >, Stokhos::MonomialProjGramSchmidtPCEBasis2< ordinal_type, value_type >, Stokhos::ProductLanczosGramSchmidtPCEBasis< ordinal_type, value_type >, Stokhos::TensorProductBasis< ordinal_type, value_type, coeff_compare_type >, Stokhos::TotalOrderBasis< ordinal_type, value_type, coeff_compare_type >, and Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >.
Referenced by Stokhos::QuadraturePseudoSpectralOperator< ordinal_t, value_t, point_compare_type >::QuadraturePseudoSpectralOperator().