44 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
45 template <
typename index_set_type>
49 const index_set_type& index_set,
51 const point_compare_type& point_compare)
53 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
58 smolyak_basis_type smolyak_basis(
59 product_basis->getCoordinateBases(), index_set, duplicate_tol);
62 smolyak_operator_type smolyak_operator(smolyak_basis,
true,
true,
68 quad_points.resize(nqp);
69 quad_weights.resize(nqp);
70 quad_values.resize(nqp);
71 typedef typename smolyak_operator_type::const_set_iterator const_iterator;
73 for (const_iterator it = smolyak_operator.set_begin();
74 it != smolyak_operator.set_end(); ++it, ++i) {
75 quad_points[i] = it->first.getTerm();
76 quad_weights[i] = it->second.first;
77 quad_values[i].resize(npc);
78 product_basis->evaluateBases(quad_points[i], quad_values[i]);
82 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
90 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
98 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
106 template <
typename ordinal_type,
typename value_type,
typename po
int_compare_type>
112 os <<
"Smolyak Sparse Grid Quadrature with " << nqp <<
" points:"
113 << std::endl <<
"Weight : Points" << std::endl;
115 os << i <<
": " << quad_weights[i] <<
" : ";
116 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_points[i].size());
118 os << quad_points[i][
j] <<
" ";
121 os <<
"Basis values at quadrature points:" << std::endl;
123 os << i <<
" " <<
": ";
124 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_values[i].size());
126 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.