Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_AnisoSparseGridQuadratureImp.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Stokhos Package
4 //
5 // Copyright 2009 NTESS and the Stokhos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "sandia_sgmga.hpp"
11 #include "sandia_rules.hpp"
12 #include "Teuchos_TimeMonitor.hpp"
13 
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,
18  ordinal_type sparse_grid_level, value_type dim_weights[],
19  value_type duplicate_tol,
20  ordinal_type growth_rate) :
21  coordinate_bases(product_basis->getCoordinateBases())
22 {
23 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
24  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::AnisoSparseGridQuadrature -- Quad Grid Generation");
25 #endif
26  ordinal_type d = product_basis->dimension();
27  ordinal_type p = product_basis->order();
28  ordinal_type sz = product_basis->size();
29  ordinal_type level = sparse_grid_level;
30 
31  // Mike's heuristic formula for computing the level
32  if (level == 0) {
33  level = static_cast<ordinal_type>(std::ceil(0.5*(p+d-1)));
34  if (level < d)
35  level = p;
36  }
37 
38  std::cout << "Sparse grid level = " << level << std::endl;
39 
40  // Compute quad points, weights, values
42 
45  for (ordinal_type i=0; i<d; i++) {
46  compute1DPoints[i] = &(getMyPoints);
47  compute1DWeights[i] = &(getMyWeights);
48  growth_rules[i] = coordinate_bases[i]->getSparseGridGrowthRule();
49  }
50 
51  // Set the static sparse grid quadrature pointer to this
52  // (this will cause a conflict if another sparse grid quadrature object
53  // is trying to access the VPISparseGrid library, but that's all we can
54  // do at this point).
55  sgq = this;
56 
57  int num_total_pts =
58  webbur::sandia_sgmga_size_total(d,&dim_weights[0], level, growth_rate,
59  &growth_rules[0]);
60 
61  ordinal_type ntot =
62  webbur::sandia_sgmga_size(d,&dim_weights[0],level,
63  &compute1DPoints[0], duplicate_tol, growth_rate,
64  &growth_rules[0]);
65 
66  Teuchos::Array<int> sparse_order(ntot*d);
67  Teuchos::Array<int> sparse_index(ntot*d);
68  Teuchos::Array<int> sparse_unique_index(num_total_pts);
69  quad_points.resize(ntot);
70  quad_weights.resize(ntot);
71  quad_values.resize(ntot);
72  Teuchos::Array<value_type> gp(ntot*d);
73 
74  webbur::sandia_sgmga_unique_index(d, &dim_weights[0], level,
75  &compute1DPoints[0],
76  duplicate_tol, ntot, num_total_pts,
77  growth_rate, &growth_rules[0],
78  &sparse_unique_index[0] );
79 
80 
81  webbur::sandia_sgmga_index(d, &dim_weights[0], level,
82  ntot, num_total_pts,
83  &sparse_unique_index[0],
84  growth_rate, &growth_rules[0],
85  &sparse_order[0], &sparse_index[0]);
86 
87  webbur::sandia_sgmga_weight(d,&dim_weights[0],level,
88  &compute1DWeights[0],
89  ntot, num_total_pts, &sparse_unique_index[0],
90  growth_rate, &growth_rules[0],
91  &quad_weights[0]);
92 
93  webbur::sandia_sgmga_point(d, &dim_weights[0], level,
94  &compute1DPoints[0],
95  ntot, &sparse_order[0], &sparse_index[0],
96  growth_rate, &growth_rules[0],
97  &gp[0]);
98 
99  for (ordinal_type i=0; i<ntot; i++) {
100  quad_values[i].resize(sz);
101  quad_points[i].resize(d);
102  for (ordinal_type j=0; j<d; j++)
103  quad_points[i][j] = gp[i*d+j];
104  product_basis->evaluateBases(quad_points[i], quad_values[i]);
105  }
106 
107  std::cout << "Number of quadrature points = " << ntot << std::endl;
108 }
109 
110 template <typename ordinal_type, typename value_type>
112 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
113 getQuadPoints() const
114 {
115  return quad_points;
116 }
117 
118 template <typename ordinal_type, typename value_type>
120 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
121 getQuadWeights() const
122 {
123  return quad_weights;
124 }
125 
126 template <typename ordinal_type, typename value_type>
128 Stokhos::AnisoSparseGridQuadrature<ordinal_type, value_type>::
129 getBasisAtQuadPoints() const
130 {
131  return quad_values;
132 }
133 
134 template <typename ordinal_type, typename value_type>
135 void
136 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
137 getMyPoints( int order, int dim, double x[] )
138 {
139  Teuchos::Array<double> quad_points;
140  Teuchos::Array<double> quad_weights;
142  ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
143  sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
144  quad_weights, quad_values);
145  TEUCHOS_ASSERT(quad_points.size() == order);
146  for (int i = 0; i<quad_points.size(); i++) {
147  x[i] = quad_points[i];
148  }
149 }
150 
151 template <typename ordinal_type, typename value_type>
152 void
153 Stokhos::AnisoSparseGridQuadrature<ordinal_type,value_type>::
154 getMyWeights( int order, int dim, double w[] )
155 {
156  Teuchos::Array<double> quad_points;
157  Teuchos::Array<double> quad_weights;
159  ordinal_type p = sgq->coordinate_bases[dim]->quadDegreeOfExactness(order);
160  sgq->coordinate_bases[dim]->getQuadPoints(p, quad_points,
161  quad_weights, quad_values);
162  TEUCHOS_ASSERT(quad_points.size() == order);
163  for (int i = 0; i<quad_points.size(); i++) {
164  w[i] = quad_weights[i];
165  }
166 }
#define TEUCHOS_FUNC_TIME_MONITOR(FUNCNAME)
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
size_type size() const
#define TEUCHOS_ASSERT(assertion_test)