12 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
13 template <
typename index_set_type>
17 const index_set_type& index_set,
19 const point_compare_type& point_compare)
21 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
26 smolyak_basis_type smolyak_basis(
27 product_basis->getCoordinateBases(), index_set, duplicate_tol);
30 smolyak_operator_type smolyak_operator(smolyak_basis,
true,
true,
36 quad_points.resize(nqp);
37 quad_weights.resize(nqp);
38 quad_values.resize(nqp);
39 typedef typename smolyak_operator_type::const_set_iterator const_iterator;
41 for (const_iterator it = smolyak_operator.set_begin();
42 it != smolyak_operator.set_end(); ++it, ++i) {
43 quad_points[i] = it->first.getTerm();
44 quad_weights[i] = it->second.first;
45 quad_values[i].resize(npc);
46 product_basis->evaluateBases(quad_points[i], quad_values[i]);
50 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
58 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
66 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
74 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
80 os <<
"Smolyak Sparse Grid Quadrature with " << nqp <<
" points:"
81 << std::endl <<
"Weight : Points" << std::endl;
83 os << i <<
": " << quad_weights[i] <<
" : ";
84 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_points[i].size());
86 os << quad_points[i][
j] <<
" ";
89 os <<
"Basis values at quadrature points:" << std::endl;
91 os << i <<
" " <<
": ";
92 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_values[i].size());
94 os << quad_values[i][
j] <<
" ";
SmolyakSparseGridQuadrature(const Teuchos::RCP< const ProductBasis< ordinal_type, value_type > > &product_basis, const index_set_type &index_set, const value_type duplicate_tol=1.0e-12, const point_compare_type &point_compare=point_compare_type())
Constructor.
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
virtual const Teuchos::Array< value_type > & getQuadWeights() const
Get quadrature weights.
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const
Get quadrature points.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const
Get values of basis at quadrature points.
virtual std::ostream & print(std::ostream &os) const
Print quadrature data.