45 template <
typename ordinal_type,
typename value_type>
49 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
64 coordinate_bases[i]->getQuadPoints(2*(coordinate_bases[i]->order()),
70 quad_weights.resize(ntot);
71 quad_values.resize(ntot);
77 quad_points[cnt].resize(d);
79 quad_values[cnt].resize(sz);
81 quad_points[cnt][
j] = gp[
j][index[
j]];
82 quad_weights[cnt] *= gw[
j][index[
j]];
88 quad_values[cnt][k] *= gv[
j][index[
j]][term[
j]];
92 while (i < d-1 && index[i] == n[i]) {
103 template <
typename ordinal_type,
typename value_type>
107 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
122 coordinate_bases[i]->getQuadPoints(quad_order,
123 gp[i], gw[i], gv[i]);
128 quad_weights.resize(ntot);
129 quad_values.resize(ntot);
135 quad_points[cnt].resize(d);
137 quad_values[cnt].resize(sz);
139 quad_points[cnt][
j] = gp[
j][index[
j]];
140 quad_weights[cnt] *= gw[
j][index[
j]];
146 quad_values[cnt][k] *= gv[
j][index[
j]][term[
j]];
150 while (i < d-1 && index[i] == n[i]) {
161 template <
typename ordinal_type,
typename value_type>
169 template <
typename ordinal_type,
typename value_type>
177 template <
typename ordinal_type,
typename value_type>
185 template <
typename ordinal_type,
typename value_type>
191 os <<
"Tensor Product Quadrature with " << nqp <<
" points:"
192 << std::endl <<
"Weight : Points" << std::endl;
194 os << i <<
": " << quad_weights[i] <<
" : ";
195 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_points[i].size());
197 os << quad_points[i][
j] <<
" ";
200 os <<
"Basis values at quadrature points:" << std::endl;
202 os << i <<
" " <<
": ";
203 for (
ordinal_type j=0; j<static_cast<ordinal_type>(quad_values[i].size());
205 os << quad_values[i][
j] <<
" ";
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
SparseArrayIterator< index_iterator, value_iterator >::value_type index(const SparseArrayIterator< index_iterator, value_iterator > &it)
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getQuadPoints() const
Get quadrature points.
virtual const Teuchos::Array< value_type > & getQuadWeights() const
Get quadrature weights.
TensorProductQuadrature(const Teuchos::RCP< const ProductBasis< ordinal_type, value_type > > &product_basis)
Constructor.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const
Get values of basis at quadrature points.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
void resize(size_type new_size, const value_type &x=value_type())
virtual std::ostream & print(std::ostream &os) const
Print quadrature data.