44 #include "sandia_sgmga.hpp"
45 #include "sandia_rules.hpp"
48 template <
typename ordinal_type,
typename value_type>
49 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
50 AnisoSparseGridQuadrature(
51 const Teuchos::RCP<
const ProductBasis<ordinal_type,value_type> >& product_basis,
55 coordinate_bases(product_basis->getCoordinateBases())
57 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
72 std::cout <<
"Sparse grid level = " << level << std::endl;
80 compute1DPoints[i] = &(getMyPoints);
81 compute1DWeights[i] = &(getMyWeights);
82 growth_rules[i] = coordinate_bases[i]->getSparseGridGrowthRule();
92 webbur::sandia_sgmga_size_total(d,&dim_weights[0], level, growth_rate,
96 webbur::sandia_sgmga_size(d,&dim_weights[0],level,
97 &compute1DPoints[0], duplicate_tol, growth_rate,
103 quad_points.resize(ntot);
104 quad_weights.resize(ntot);
105 quad_values.resize(ntot);
108 webbur::sandia_sgmga_unique_index(d, &dim_weights[0], level,
110 duplicate_tol, ntot, num_total_pts,
111 growth_rate, &growth_rules[0],
112 &sparse_unique_index[0] );
115 webbur::sandia_sgmga_index(d, &dim_weights[0], level,
117 &sparse_unique_index[0],
118 growth_rate, &growth_rules[0],
119 &sparse_order[0], &sparse_index[0]);
121 webbur::sandia_sgmga_weight(d,&dim_weights[0],level,
122 &compute1DWeights[0],
123 ntot, num_total_pts, &sparse_unique_index[0],
124 growth_rate, &growth_rules[0],
127 webbur::sandia_sgmga_point(d, &dim_weights[0], level,
129 ntot, &sparse_order[0], &sparse_index[0],
130 growth_rate, &growth_rules[0],
134 quad_values[i].resize(sz);
135 quad_points[i].resize(d);
137 quad_points[i][
j] = gp[i*d+
j];
138 product_basis->evaluateBases(quad_points[i], quad_values[i]);
141 std::cout <<
"Number of quadrature points = " << ntot << std::endl;
144 template <
typename ordinal_type,
typename value_type>
146 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
147 getQuadPoints()
const
152 template <
typename ordinal_type,
typename value_type>
154 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
155 getQuadWeights()
const
160 template <
typename ordinal_type,
typename value_type>
162 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
163 getBasisAtQuadPoints()
const
168 template <
typename ordinal_type,
typename value_type>
170 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
171 getMyPoints(
int order,
int dim,
double x[] )
176 ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
177 sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
178 quad_weights, quad_values);
180 for (
int i = 0; i<quad_points.
size(); i++) {
181 x[i] = quad_points[i];
185 template <
typename ordinal_type,
typename value_type>
187 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
188 getMyWeights(
int order,
int dim,
double w[] )
193 ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
194 sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
195 quad_weights, quad_values);
197 for (
int i = 0; i<quad_points.
size(); i++) {
198 w[i] = quad_weights[i];
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
ordinal_type sz
Expansions size.
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
#define TEUCHOS_ASSERT(assertion_test)