| Stokhos Package Browser (Single Doxygen Collection)
    Version of the Day
    | 
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.  More... | |
| ordinal_type | getNumSmolyakTerms () const | 
| Return number of terms in Smolyak formula.  More... | |
| Teuchos::RCP< const tensor_product_basis_type > | getTensorProductBasis (ordinal_type i) const | 
| Return ith tensor product basis.  More... | |
| ordinal_type | getSmolyakCoefficient (ordinal_type i) const | 
| Return ith smolyak coefficient.  More... | |
| 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) | |
|  Public Member Functions inherited from Stokhos::ProductBasis< ordinal_type, value_type > | |
| ProductBasis () | |
| Constructor.  More... | |
| virtual | ~ProductBasis () | 
| Destructor.  More... | |
|  Public Member Functions inherited from Stokhos::OrthogPolyBasis< ordinal_type, value_type > | |
| OrthogPolyBasis () | |
| Constructor.  More... | |
| virtual | ~OrthogPolyBasis () | 
| Destructor.  More... | |
| 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.  More... | |
| ordinal_type | p | 
| Total order of basis.  More... | |
| ordinal_type | d | 
| Total dimension of basis.  More... | |
| ordinal_type | sz | 
| Total size of basis.  More... | |
| Teuchos::Array< Teuchos::RCP < const OneDOrthogPolyBasis < ordinal_type, value_type > > > | bases | 
| Array of bases.  More... | |
| value_type | sparse_tol | 
| Tolerance for computing sparse Cijk.  More... | |
| coeff_type | max_orders | 
| Maximum orders for each dimension.  More... | |
| coeff_set_type | basis_set | 
| Basis set.  More... | |
| coeff_map_type | basis_map | 
| Basis map.  More... | |
| Teuchos::Array< ordinal_type > | smolyak_coeffs | 
| Smolyak coefficients.  More... | |
| Teuchos::Array< value_type > | norms | 
| Norms.  More... | |
| Teuchos::Array< Teuchos::Array < value_type > > | basis_eval_tmp | 
| Temporary array used in basis evaluation.  More... | |
| Teuchos::Array< Teuchos::RCP < tensor_product_basis_type > > | tp_bases | 
| Tensor product bases comprising Smolyak set.  More... | |
| SmolyakPredicate < TensorProductPredicate < ordinal_type > > | sm_pred | 
| Predicate for building sparse triple products.  More... | |
| Private Member Functions | |
| SmolyakBasis (const SmolyakBasis &) | |
| SmolyakBasis & | operator= (const SmolyakBasis &b) | 
| Implementation of Stokhos::OrthogPolyBasis methods | |
| ordinal_type | order () const | 
| Return order of basis.  More... | |
| ordinal_type | dimension () const | 
| Return dimension of basis.  More... | |
| virtual ordinal_type | size () const | 
| Return total size of basis.  More... | |
| 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.  More... | |
| 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.  More... | |
| virtual value_type | evaluateZero (ordinal_type i) const | 
| Evaluate basis polynomial iat zero.  More... | |
| 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.  More... | |
| virtual const std::string & | getName () const | 
| Return string name of basis.  More... | |
| 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.  More... | |
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
Definition at line 30 of file Stokhos_SmolyakBasis.hpp.
| typedef MultiIndex<ordinal_type> Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >::coeff_type | 
Definition at line 140 of file Stokhos_SmolyakBasis.hpp.
| typedef TensorProductBasis<ordinal_type,value_type, LexographicLess<coeff_type> > Stokhos::SmolyakBasis< ordinal_type, value_type, coeff_compare_type >::tensor_product_basis_type | 
Definition at line 143 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Definition at line 167 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Definition at line 168 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Definition at line 169 of file Stokhos_SmolyakBasis.hpp.
| 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 | 
Destructor.
Definition at line 181 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | private | 
| 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_, | ||
| const ordering_type & | coeff_compare | ||
| ) | 
Definition at line 16 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return order of basis.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 188 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return dimension of basis.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 196 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return total size of basis.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 204 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return array storing norm-squared of each basis polynomial.
Entry  of returned array is given by
 of returned array is given by  for
 for  where
 where  is size()-1.
 is size()-1. 
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 212 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return norm squared of basis polynomial i. 
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 220 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Compute triple product tensor.
The  entry of the tensor
 entry of the tensor  is given by
 is given by  where
 where  represents basis polynomial
 represents basis polynomial  and
 and  where
 where  is size()-1.
 is size()-1. 
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 228 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Compute linear triple product tensor where k = 0,1,..,d.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 241 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Evaluate basis polynomial i at zero. 
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 260 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | 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 >.
Definition at line 274 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Print basis to stream os. 
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 292 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return string name of basis.
Implements Stokhos::OrthogPolyBasis< ordinal_type, value_type >.
Definition at line 326 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Get orders of each coordinate polynomial given an index i. 
The returned array is of size  , where
, where  is the dimension of the basis, and entry
 is the dimension of the basis, and entry  is given by
 is given by  where
 where  .
. 
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Definition at line 307 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Get index of the multivariate polynomial given orders of each coordinate.
Given the array term storing  , returns the index
, returns the index  such that
 such that  .
. 
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Definition at line 315 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return coordinate bases.
Array is of size dimension().
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Definition at line 334 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | virtual | 
Return maximum order allowable for each coordinate basis.
Implements Stokhos::ProductBasis< ordinal_type, value_type >.
Definition at line 342 of file Stokhos_SmolyakBasisImp.hpp.
| 
 | inline | 
Return number of terms in Smolyak formula.
Definition at line 146 of file Stokhos_SmolyakBasis.hpp.
| 
 | inline | 
Return ith tensor product basis.
Definition at line 150 of file Stokhos_SmolyakBasis.hpp.
| 
 | inline | 
Return ith smolyak coefficient.
Definition at line 153 of file Stokhos_SmolyakBasis.hpp.
| 
 | private | 
| 
 | protected | 
Name of basis.
Definition at line 172 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Total order of basis.
Definition at line 175 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Total dimension of basis.
Definition at line 178 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Total size of basis.
Definition at line 181 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Array of bases.
Definition at line 184 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Tolerance for computing sparse Cijk.
Definition at line 187 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Maximum orders for each dimension.
Definition at line 190 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Basis set.
Definition at line 193 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Basis map.
Definition at line 196 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Smolyak coefficients.
Definition at line 199 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Norms.
Definition at line 202 of file Stokhos_SmolyakBasis.hpp.
| 
 | mutableprotected | 
Temporary array used in basis evaluation.
Definition at line 205 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Tensor product bases comprising Smolyak set.
Definition at line 208 of file Stokhos_SmolyakBasis.hpp.
| 
 | protected | 
Predicate for building sparse triple products.
Definition at line 226 of file Stokhos_SmolyakBasis.hpp.
 1.8.5
 1.8.5