10 #include "sandia_sgmga.hpp"
11 #include "sandia_rules.hpp"
14 template <
typename ordinal_type,
typename value_type>
15 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
16 AnisoSparseGridQuadrature(
17 const Teuchos::RCP<
const ProductBasis<ordinal_type,value_type> >& product_basis,
21 coordinate_bases(product_basis->getCoordinateBases())
23 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
38 std::cout <<
"Sparse grid level = " << level << std::endl;
46 compute1DPoints[i] = &(getMyPoints);
47 compute1DWeights[i] = &(getMyWeights);
48 growth_rules[i] = coordinate_bases[i]->getSparseGridGrowthRule();
58 webbur::sandia_sgmga_size_total(d,&dim_weights[0], level, growth_rate,
62 webbur::sandia_sgmga_size(d,&dim_weights[0],level,
63 &compute1DPoints[0], duplicate_tol, growth_rate,
69 quad_points.resize(ntot);
70 quad_weights.resize(ntot);
71 quad_values.resize(ntot);
74 webbur::sandia_sgmga_unique_index(d, &dim_weights[0], level,
76 duplicate_tol, ntot, num_total_pts,
77 growth_rate, &growth_rules[0],
78 &sparse_unique_index[0] );
81 webbur::sandia_sgmga_index(d, &dim_weights[0], level,
83 &sparse_unique_index[0],
84 growth_rate, &growth_rules[0],
85 &sparse_order[0], &sparse_index[0]);
87 webbur::sandia_sgmga_weight(d,&dim_weights[0],level,
89 ntot, num_total_pts, &sparse_unique_index[0],
90 growth_rate, &growth_rules[0],
93 webbur::sandia_sgmga_point(d, &dim_weights[0], level,
95 ntot, &sparse_order[0], &sparse_index[0],
96 growth_rate, &growth_rules[0],
100 quad_values[i].resize(sz);
101 quad_points[i].resize(d);
103 quad_points[i][
j] = gp[i*d+
j];
104 product_basis->evaluateBases(quad_points[i], quad_values[i]);
107 std::cout <<
"Number of quadrature points = " << ntot << std::endl;
110 template <
typename ordinal_type,
typename value_type>
112 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
113 getQuadPoints()
const
118 template <
typename ordinal_type,
typename value_type>
120 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
121 getQuadWeights()
const
126 template <
typename ordinal_type,
typename value_type>
128 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
129 getBasisAtQuadPoints()
const
134 template <
typename ordinal_type,
typename value_type>
136 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
137 getMyPoints(
int order,
int dim,
double x[] )
142 ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
143 sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
144 quad_weights, quad_values);
146 for (
int i = 0; i<quad_points.
size(); i++) {
147 x[i] = quad_points[i];
151 template <
typename ordinal_type,
typename value_type>
153 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
154 getMyWeights(
int order,
int dim,
double w[] )
159 ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
160 sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
161 quad_weights, quad_values);
163 for (
int i = 0; i<quad_points.
size(); i++) {
164 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)