42#ifndef STOKHOS_SMOLYAK_BASIS_HPP
43#define STOKHOS_SMOLYAK_BASIS_HPP
47#include "Teuchos_RCP.hpp"
60 typename coeff_compare_type =
61 TotalOrderLess<MultiIndex<ordinal_type> > >
72 template <
typename index_set_type>
75 value_type> > >&
bases,
76 const index_set_type& index_set,
78 const coeff_compare_type& coeff_compare = coeff_compare_type());
87 ordinal_type
order()
const;
93 virtual ordinal_type
size()
const;
100 virtual const Teuchos::Array<value_type>&
norm_squared()
const;
103 virtual const value_type&
norm_squared(ordinal_type i)
const;
113 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
118 Teuchos::RCP< Stokhos::Sparse3Tensor<ordinal_type, value_type> >
130 const Teuchos::ArrayView<const value_type>& point,
131 Teuchos::Array<value_type>& basis_vals)
const;
134 virtual void print(std::ostream& os)
const;
137 virtual const std::string&
getName()
const;
150 virtual const MultiIndex<ordinal_type>&
term(ordinal_type i)
const;
157 virtual ordinal_type
index(
const MultiIndex<ordinal_type>&
term)
const;
174 LexographicLess<coeff_type>
181 Teuchos::RCP<const tensor_product_basis_type>
216 Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type, value_type> > >
bases;
240 Teuchos::Array< Teuchos::RCP< tensor_product_basis_type > >
tp_bases;
243 template <
typename tp_predicate_type>
249 for (ordinal_type i=0; i<
tp_preds.size(); ++i)
258 SmolyakPredicate< TensorProductPredicate<ordinal_type> >
sm_pred;
Abstract base class for 1-D orthogonal polynomials.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
Teuchos::Array< Teuchos::RCP< tensor_product_basis_type > > tp_bases
Tensor product bases comprising Smolyak set.
virtual void evaluateBases(const Teuchos::ArrayView< const value_type > &point, Teuchos::Array< value_type > &basis_vals) const
Evaluate basis polynomials at given point point.
virtual const Teuchos::Array< value_type > & norm_squared() const
Return array storing norm-squared of each basis polynomial.
ordinal_type dimension() const
Return dimension of basis.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeLinearTripleProductTensor() const
Compute linear triple product tensor where k = 0,1,..,d.
Teuchos::Array< coeff_type > coeff_map_type
virtual ordinal_type index(const MultiIndex< ordinal_type > &term) const
Get index of the multivariate polynomial given orders of each coordinate.
Teuchos::Array< Teuchos::Array< value_type > > basis_eval_tmp
Temporary array used in basis evaluation.
Teuchos::Array< ordinal_type > smolyak_coeffs
Smolyak coefficients.
TensorProductBasis< ordinal_type, value_type, LexographicLess< coeff_type > > tensor_product_basis_type
std::map< coeff_type, ordinal_type, coeff_compare_type > coeff_set_type
std::string name
Name of basis.
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
ordinal_type d
Total dimension of basis.
ordinal_type order() const
Return order of basis.
coeff_map_type basis_map
Basis map.
MultiIndex< ordinal_type > coeff_type
virtual ordinal_type size() const
Return total size of basis.
coeff_type max_orders
Maximum orders for each dimension.
Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > bases
Array of bases.
ordinal_type sz
Total size of basis.
value_type sparse_tol
Tolerance for computing sparse Cijk.
SmolyakBasis(const SmolyakBasis &)
Teuchos::RCP< const tensor_product_basis_type > getTensorProductBasis(ordinal_type i) const
Return ith tensor product basis.
coeff_set_type basis_set
Basis set.
virtual const std::string & getName() const
Return string name of basis.
Teuchos::Array< value_type > norms
Norms.
virtual MultiIndex< ordinal_type > getMaxOrders() const
Return maximum order allowable for each coordinate basis.
SmolyakPredicate< TensorProductPredicate< ordinal_type > > sm_pred
Predicate for building sparse triple products.
ordinal_type getNumSmolyakTerms() const
Return number of terms in Smolyak formula.
MultiIndex< ordinal_type > multiindex_type
Teuchos::Array< Teuchos::RCP< const OneDOrthogPolyBasis< ordinal_type, value_type > > > getCoordinateBases() const
Return coordinate bases.
SmolyakBasis & operator=(const SmolyakBasis &b)
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.
virtual const MultiIndex< ordinal_type > & term(ordinal_type i) const
Get orders of each coordinate polynomial given an index i.
virtual void print(std::ostream &os) const
Print basis to stream os.
virtual ~SmolyakBasis()
Destructor.
ordinal_type p
Total order of basis.
ordinal_type getSmolyakCoefficient(ordinal_type i) const
Return ith smolyak coefficient.
virtual value_type evaluateZero(ordinal_type i) const
Evaluate basis polynomial i at zero.
Multivariate orthogonal polynomial basis generated from a tensor product of univariate polynomials.
Top-level namespace for Stokhos classes and functions.
Predicate functor for building sparse triple products.
bool operator()(const coeff_type &term) const
Teuchos::Array< tp_predicate_type > tp_preds