|
Stokhos
Development
|
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid. More...
#include <Stokhos_SmolyakBasis.hpp>


Classes | |
| struct | SmolyakPredicate |
| Predicate functor for building sparse triple products. More... | |
Public Types | |
| typedef MultiIndex< ordinal_type > | coeff_type |
|
typedef TensorProductBasis < ordinal_type, value_type, LexographicLess< coeff_type > > | tensor_product_basis_type |
Public Member Functions | |
| template<typename index_set_type > | |
| SmolyakBasis (const Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > &bases, const index_set_type &index_set, const value_type &sparse_tol=1.0e-12, const coeff_compare_type &coeff_compare=coeff_compare_type()) | |
| Constructor. More... | |
| virtual | ~SmolyakBasis () |
| Destructor. | |
| ordinal_type | getNumSmolyakTerms () const |
| Return number of terms in Smolyak formula. | |
|
Teuchos::RCP< const tensor_product_basis_type > | getTensorProductBasis (ordinal_type i) const |
| Return ith tensor product basis. | |
| ordinal_type | getSmolyakCoefficient (ordinal_type i) const |
| Return ith smolyak coefficient. | |
| template<typename index_set_type > | |
| SmolyakBasis (const Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > &bases_, const index_set_type &index_set, const value_type &sparse_tol_, const ordering_type &coeff_compare) | |
Implementation of Stokhos::OrthogPolyBasis methods | |
| ordinal_type | order () const |
| Return order of basis. | |
| ordinal_type | dimension () const |
| Return dimension of basis. | |
| virtual ordinal_type | size () const |
| Return total size of basis. | |
| virtual const Teuchos::Array < value_type > & | norm_squared () const |
| Return array storing norm-squared of each basis polynomial. More... | |
| virtual const value_type & | norm_squared (ordinal_type i) const |
Return norm squared of basis polynomial i. | |
| virtual Teuchos::RCP < Stokhos::Sparse3Tensor < ordinal_type, value_type > > | computeTripleProductTensor () const |
| Compute triple product tensor. More... | |
|
virtual Teuchos::RCP < Stokhos::Sparse3Tensor < ordinal_type, value_type > > | computeLinearTripleProductTensor () const |
| Compute linear triple product tensor where k = 0,1,..,d. | |
| virtual value_type | evaluateZero (ordinal_type i) const |
Evaluate basis polynomial i at zero. | |
| virtual void | evaluateBases (const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const |
Evaluate basis polynomials at given point point. More... | |
| virtual void | print (std::ostream &os) const |
Print basis to stream os. | |
| virtual const std::string & | getName () const |
| Return string name of basis. | |
Implementation of Stokhos::ProductBasis methods | |
| virtual const MultiIndex < ordinal_type > & | term (ordinal_type i) const |
Get orders of each coordinate polynomial given an index i. More... | |
| virtual ordinal_type | index (const MultiIndex< ordinal_type > &term) const |
| Get index of the multivariate polynomial given orders of each coordinate. More... | |
| Teuchos::Array< Teuchos::RCP < const OneDOrthogPolyBasis < ordinal_type, value_type > > > | getCoordinateBases () const |
| Return coordinate bases. More... | |
| virtual MultiIndex< ordinal_type > | getMaxOrders () const |
| Return maximum order allowable for each coordinate basis. | |
Public Member Functions inherited from Stokhos::ProductBasis< ordinal_type, value_type > | |
| ProductBasis () | |
| Constructor. | |
| virtual | ~ProductBasis () |
| Destructor. | |
Public Member Functions inherited from Stokhos::OrthogPolyBasis< ordinal_type, value_type > | |
| OrthogPolyBasis () | |
| Constructor. | |
| virtual | ~OrthogPolyBasis () |
| Destructor. | |
Protected Types | |
|
typedef std::map< coeff_type, ordinal_type, coeff_compare_type > | coeff_set_type |
|
typedef Teuchos::Array < coeff_type > | coeff_map_type |
| typedef MultiIndex< ordinal_type > | multiindex_type |
Protected Attributes | |
| std::string | name |
| Name of basis. | |
| ordinal_type | p |
| Total order of basis. | |
| ordinal_type | d |
| Total dimension of basis. | |
| ordinal_type | sz |
| Total size of basis. | |
|
Teuchos::Array< Teuchos::RCP < const OneDOrthogPolyBasis < ordinal_type, value_type > > > | bases |
| Array of bases. | |
| value_type | sparse_tol |
| Tolerance for computing sparse Cijk. | |
| coeff_type | max_orders |
| Maximum orders for each dimension. | |
| coeff_set_type | basis_set |
| Basis set. | |
| coeff_map_type | basis_map |
| Basis map. | |
| Teuchos::Array< ordinal_type > | smolyak_coeffs |
| Smolyak coefficients. | |
| Teuchos::Array< value_type > | norms |
| Norms. | |
|
Teuchos::Array< Teuchos::Array < value_type > > | basis_eval_tmp |
| Temporary array used in basis evaluation. | |
|
Teuchos::Array< Teuchos::RCP < tensor_product_basis_type > > | tp_bases |
| Tensor product bases comprising Smolyak set. | |
|
SmolyakPredicate < TensorProductPredicate < ordinal_type > > | sm_pred |
| Predicate for building sparse triple products. | |
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
| Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >::SmolyakBasis | ( | const Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > & | bases, |
| const index_set_type & | index_set, | ||
| const value_type & | sparse_tol = 1.0e-12, |
||
| const coeff_compare_type & | coeff_compare = coeff_compare_type() |
||
| ) |
Constructor.
| bases | array of 1-D coordinate bases |
| sparse_tol | tolerance used to drop terms in sparse triple-product tensors |
|
virtual |
Compute triple product tensor.
The
entry of the tensor
is given by
where
represents basis polynomial
and
where
is size()-1.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
|
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.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Referenced by Stokhos::SmolyakPseudoSpectralOperator< ordinal_t, value_t, point_compare_type >::SmolyakPseudoSpectralOperator().
|
virtual |
Return coordinate bases.
Array is of size dimension().
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Referenced by Stokhos::SmolyakPseudoSpectralOperator< ordinal_t, value_t, point_compare_type >::SmolyakPseudoSpectralOperator().
|
virtual |
Get index of the multivariate polynomial given orders of each coordinate.
Given the array term storing
, returns the index
such that
.
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Referenced by Stokhos::SmolyakPseudoSpectralOperator< ordinal_t, value_t, point_compare_type >::SmolyakPseudoSpectralOperator().
|
virtual |
Return array storing norm-squared of each basis polynomial.
Entry
of returned array is given by
for
where
is size()-1.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Referenced by Stokhos::SmolyakPseudoSpectralOperator< ordinal_t, value_t, point_compare_type >::SmolyakPseudoSpectralOperator().
|
virtual |
Get orders of each coordinate polynomial given an index i.
The returned array is of size
, where
is the dimension of the basis, and entry
is given by
where
.
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
1.8.5