Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_TensorProductQuadratureImp.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 "Teuchos_TimeMonitor.hpp"
11 
12 template <typename ordinal_type, typename value_type>
15 {
16 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
17  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::TensorProductQuadrature -- Quad Grid Generation");
18 #endif
19  ordinal_type d = product_basis->dimension();
20  ordinal_type sz = product_basis->size();
21 
22  Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type,value_type> > > coordinate_bases = product_basis->getCoordinateBases();
23 
24  // Compute quad points, weights, values
29  ordinal_type ntot = 1;
30  for (ordinal_type i=0; i<d; i++) {
31  coordinate_bases[i]->getQuadPoints(2*(coordinate_bases[i]->order()),
32  gp[i], gw[i], gv[i]);
33  n[i] = gp[i].size();
34  ntot *= n[i];
35  }
36  quad_points.resize(ntot);
37  quad_weights.resize(ntot);
38  quad_values.resize(ntot);
40  for (ordinal_type i=0; i<d; i++)
41  index[i] = 0;
42  ordinal_type cnt = 0;
43  while (cnt < ntot) {
44  quad_points[cnt].resize(d);
45  quad_weights[cnt] = value_type(1.0);
46  quad_values[cnt].resize(sz);
47  for (ordinal_type j=0; j<d; j++) {
48  quad_points[cnt][j] = gp[j][index[j]];
49  quad_weights[cnt] *= gw[j][index[j]];
50  }
51  for (ordinal_type k=0; k<sz; k++) {
52  quad_values[cnt][k] = value_type(1.0);
53  const MultiIndex<ordinal_type>& term = product_basis->term(k);
54  for (ordinal_type j=0; j<d; j++)
55  quad_values[cnt][k] *= gv[j][index[j]][term[j]];
56  }
57  ++index[0];
58  ordinal_type i = 0;
59  while (i < d-1 && index[i] == n[i]) {
60  index[i] = 0;
61  ++i;
62  ++index[i];
63  }
64  ++cnt;
65  }
66 
67  //std::cout << "Number of quadrature points = " << ntot << std::endl;
68 }
69 
70 template <typename ordinal_type, typename value_type>
73 {
74 #ifdef STOKHOS_TEUCHOS_TIME_MONITOR
75  TEUCHOS_FUNC_TIME_MONITOR("Stokhos::TensorProductQuadrature -- Quad Grid Generation");
76 #endif
77  ordinal_type d = product_basis->dimension();
78  ordinal_type sz = product_basis->size();
79 
80  Teuchos::Array< Teuchos::RCP<const OneDOrthogPolyBasis<ordinal_type,value_type> > > coordinate_bases = product_basis->getCoordinateBases();
81 
82  // Compute quad points, weights, values
87  ordinal_type ntot = 1;
88  for (ordinal_type i=0; i<d; i++) {
89  coordinate_bases[i]->getQuadPoints(quad_order,
90  gp[i], gw[i], gv[i]);
91  n[i] = gp[i].size();
92  ntot *= n[i];
93  }
94  quad_points.resize(ntot);
95  quad_weights.resize(ntot);
96  quad_values.resize(ntot);
98  for (ordinal_type i=0; i<d; i++)
99  index[i] = 0;
100  ordinal_type cnt = 0;
101  while (cnt < ntot) {
102  quad_points[cnt].resize(d);
103  quad_weights[cnt] = value_type(1.0);
104  quad_values[cnt].resize(sz);
105  for (ordinal_type j=0; j<d; j++) {
106  quad_points[cnt][j] = gp[j][index[j]];
107  quad_weights[cnt] *= gw[j][index[j]];
108  }
109  for (ordinal_type k=0; k<sz; k++) {
110  quad_values[cnt][k] = value_type(1.0);
111  MultiIndex<ordinal_type> term = product_basis->term(k);
112  for (ordinal_type j=0; j<d; j++)
113  quad_values[cnt][k] *= gv[j][index[j]][term[j]];
114  }
115  ++index[0];
116  ordinal_type i = 0;
117  while (i < d-1 && index[i] == n[i]) {
118  index[i] = 0;
119  ++i;
120  ++index[i];
121  }
122  ++cnt;
123  }
124 
125  //std::cout << "Number of quadrature points = " << ntot << std::endl;
126 }
127 
128 template <typename ordinal_type, typename value_type>
132 {
133  return quad_points;
134 }
135 
136 template <typename ordinal_type, typename value_type>
140 {
141  return quad_weights;
142 }
143 
144 template <typename ordinal_type, typename value_type>
148 {
149  return quad_values;
150 }
151 
152 template <typename ordinal_type, typename value_type>
153 std::ostream&
155 print(std::ostream& os) const
156 {
157  ordinal_type nqp = quad_weights.size();
158  os << "Tensor Product Quadrature with " << nqp << " points:"
159  << std::endl << "Weight : Points" << std::endl;
160  for (ordinal_type i=0; i<nqp; i++) {
161  os << i << ": " << quad_weights[i] << " : ";
162  for (ordinal_type j=0; j<static_cast<ordinal_type>(quad_points[i].size());
163  j++)
164  os << quad_points[i][j] << " ";
165  os << std::endl;
166  }
167  os << "Basis values at quadrature points:" << std::endl;
168  for (ordinal_type i=0; i<nqp; i++) {
169  os << i << " " << ": ";
170  for (ordinal_type j=0; j<static_cast<ordinal_type>(quad_values[i].size());
171  j++)
172  os << quad_values[i][j] << " ";
173  os << std::endl;
174  }
175 
176  return os;
177 }
#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.
size_type size() const
int n