43#include "Teuchos_TimeMonitor.hpp"
45template <
typename ordinal_type,
typename value_type>
49#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
50 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::TensorProductQuadrature -- Quad Grid Generation");
55 Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type,value_type> > > coordinate_bases = product_basis->getCoordinateBases();
58 Teuchos::Array< Teuchos::Array<value_type> > gp(d);
59 Teuchos::Array< Teuchos::Array<value_type> > gw(d);
60 Teuchos::Array< Teuchos::Array< Teuchos::Array<value_type> > > gv(d);
61 Teuchos::Array<ordinal_type> n(d);
64 coordinate_bases[i]->getQuadPoints(2*(coordinate_bases[i]->order()),
72 Teuchos::Array<ordinal_type>
index(d);
92 while (i < d-1 &&
index[i] == n[i]) {
103template <
typename ordinal_type,
typename value_type>
107#ifdef STOKHOS_TEUCHOS_TIME_MONITOR
108 TEUCHOS_FUNC_TIME_MONITOR(
"Stokhos::TensorProductQuadrature -- Quad Grid Generation");
113 Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type,value_type> > > coordinate_bases = product_basis->getCoordinateBases();
116 Teuchos::Array< Teuchos::Array<value_type> > gp(d);
117 Teuchos::Array< Teuchos::Array<value_type> > gw(d);
118 Teuchos::Array< Teuchos::Array< Teuchos::Array<value_type> > > gv(d);
119 Teuchos::Array<ordinal_type> n(d);
122 coordinate_bases[i]->getQuadPoints(quad_order,
123 gp[i], gw[i], gv[i]);
130 Teuchos::Array<ordinal_type>
index(d);
150 while (i < d-1 &&
index[i] == n[i]) {
161template <
typename ordinal_type,
typename value_type>
162const Teuchos::Array< Teuchos::Array<value_type> >&
169template <
typename ordinal_type,
typename value_type>
170const Teuchos::Array<value_type>&
177template <
typename ordinal_type,
typename value_type>
178const Teuchos::Array< Teuchos::Array<value_type> >&
185template <
typename ordinal_type,
typename value_type>
188print(std::ostream& os)
const
191 os <<
"Tensor Product Quadrature with " << nqp <<
" points:"
192 << std::endl <<
"Weight : Points" << std::endl;
200 os <<
"Basis values at quadrature points:" << std::endl;
202 os << i <<
" " <<
": ";
A multidimensional index.
Abstract base class for multivariate orthogonal polynomials generated from tensor products of univari...
ordinal_type size() const
Return size.
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.
Teuchos::Array< Teuchos::Array< value_type > > quad_values
Quadrature values.
virtual std::ostream & print(std::ostream &os) const
Print quadrature data.
TensorProductQuadrature(const Teuchos::RCP< const ProductBasis< ordinal_type, value_type > > &product_basis)
Constructor.
virtual ordinal_type size() const
Get number of quadrature points.
Teuchos::Array< value_type > quad_weights
Quadrature weights.
virtual const Teuchos::Array< Teuchos::Array< value_type > > & getBasisAtQuadPoints() const
Get values of basis at quadrature points.
Teuchos::Array< Teuchos::Array< value_type > > quad_points
Quadrature points.