Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp
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 
14 #include "Teuchos_ArrayView.hpp"
15 
16 #include "Stokhos.hpp"
18 
19 namespace SmolyakBasisUtilsUnitTest {
20 
21  // Common setup for unit tests
22  template <typename ordinal_type, typename value_type>
23  struct UnitTestSetup {
27 
29  rtol = 1e-12;
30  atol = 1e-12;
31  apply_rtol = 1e-2;
32  apply_atol = 1e-3;
33  p = 4;
34  d = 2;
35  }
36 
37  };
38 
39  typedef int ordinal_type;
40  typedef double value_type;
43 
44  // Function for testing quadratures
46  value_type val = 0.0;
47  for (int i=0; i<x.size(); i++)
48  val += x[i];
49  return std::exp(val);
50  }
51 
52  // Function for testing quadratures
54  value_type val = 0.0;
55  for (int i=0; i<x.size(); i++)
56  val += x[i];
57  return std::sin(val);
58  }
59 
60  TEUCHOS_UNIT_TEST( Direct, Linear ) {
61  // Build isotropic Smolyak basis of dimension d and order p with
62  // linear growth
64  for (ordinal_type i=0; i<setup.d; i++)
68 
69  // Build corresponding pseudospectral operator
71  *smolyak_basis, false);
72 
73  // Generate sparse grids using original approach
74  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
75  smolyak_basis, setup.p, 1e-12, Pecos::SLOW_RESTRICTED_GROWTH);
77  *smolyak_basis, quad);
78 
79  // Test grid
81  sm_op, quad_op, setup.rtol, setup.atol, out);
82 
83  // Test apply
84  success = success &&
86  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
87 
88  // Test transpose apply
89  success = success &&
91  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
92 
93  // Test discrete orthogonality
94  success = success &&
96  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
97  }
98 
99  TEUCHOS_UNIT_TEST( Direct, ModerateLinear ) {
100  // Build isotropic Smolyak basis of dimension d and order p with
101  // moderate linear growth
103  for (ordinal_type i=0; i<setup.d; i++)
107 
108  // Build corresponding pseudospectral operator
110  *smolyak_basis, false);
111 
112  // Generate sparse grids using original approach
113  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
114  smolyak_basis, setup.p, 1e-12, Pecos::MODERATE_RESTRICTED_GROWTH);
116  *smolyak_basis, quad);
117 
118  // Test grid
120  sm_op, quad_op, setup.rtol, setup.atol, out);
121 
122  // Test apply
123  success = success &&
125  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
126 
127  // Test transpose apply
128  success = success &&
130  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
131 
132  // Direct apply will not satisfy discrete orthogonality
133  }
134 
135  TEUCHOS_UNIT_TEST( Smolyak, Linear ) {
136  // Build isotropic Smolyak basis of dimension d and order p with
137  // linear growth
139  for (ordinal_type i=0; i<setup.d; i++)
143 
144  // Build corresponding pseudospectral operator
146  *smolyak_basis, true);
147 
148  // Generate sparse grids using original approach
149  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
150  smolyak_basis, setup.p, 1e-12, Pecos::SLOW_RESTRICTED_GROWTH);
152  *smolyak_basis, quad);
153 
155  smolyak_basis);
156  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
157 
158  // Test grid
160  sm_op, quad_op, setup.rtol, setup.atol, out);
161 
162  // Test apply
163  success = success &&
165  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
166 
167  // Test transpose apply
168  success = success &&
170  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
171 
172  // Test discrete orthogonality
173  success = success &&
175  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
176  }
177 
178  TEUCHOS_UNIT_TEST( Smolyak, ModerateLinear ) {
179  // Build isotropic Smolyak basis of dimension d and order p with
180  // moderate linear growth
182  for (ordinal_type i=0; i<setup.d; i++)
186 
187  // Build corresponding pseudospectral operator
189  *smolyak_basis, true);
190 
191  // Generate sparse grids using original approach
192  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
193  smolyak_basis, setup.p, 1e-12, Pecos::MODERATE_RESTRICTED_GROWTH);
195  *smolyak_basis, quad);
196 
198  smolyak_basis);
199  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
200 
201  // Test grid
203  sm_op, quad_op, setup.rtol, setup.atol, out);
204 
205  // Test apply
206  success = success &&
208  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
209 
210  // Test transpose apply
211  success = success &&
213  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
214 
215  // Test discrete orthogonality
216  success = success &&
218  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
219  }
220 
221 #ifdef HAVE_STOKHOS_DAKOTA
222 
223  TEUCHOS_UNIT_TEST( Direct, ClenshawCurtis ) {
224  // Build isotropic Smolyak basis of dimension d and order p with
225  // exponential growth
227  for (ordinal_type i=0; i<setup.d; i++)
231 
232  // Build corresponding pseudospectral operator
234  *smolyak_basis, false);
235 
236  // Generate sparse grids using original approach
237  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
238  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
240  *smolyak_basis, quad);
241 
242  // Test grid
244  sm_op, quad_op, setup.rtol, setup.atol, out);
245 
246  // Test apply
247  success = success &&
249  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
250 
251  // Test transpose apply
252  success = success &&
254  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
255 
256  // Direct apply will not satisfy discrete orthogonality
257  }
258 
259  TEUCHOS_UNIT_TEST( Direct, GaussPatterson ) {
260  // Build isotropic Smolyak basis of dimension d and order p with
261  // exponential Gauss-Patterson growth
263  for (ordinal_type i=0; i<setup.d; i++)
267 
268  // Build corresponding pseudospectral operator
270  *smolyak_basis, false);
271 
272  // Generate sparse grids using original approach
273  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
274  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
276  *smolyak_basis, quad);
277 
278  // Test grid
280  sm_op, quad_op, setup.rtol, setup.atol, out);
281 
282  // Test apply
283  success = success &&
285  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
286 
287  // Test transpose apply
288  success = success &&
290  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
291 
292  // Direct apply will not satisfy discrete orthogonality
293  }
294 
295  TEUCHOS_UNIT_TEST( Smolyak, ClenshawCurtis ) {
296  // Build isotropic Smolyak basis of dimension d and order p with
297  // exponential growth
299  for (ordinal_type i=0; i<setup.d; i++)
303 
304  // Build corresponding pseudospectral operator
306  *smolyak_basis, true);
307 
308  // Generate sparse grids using original approach
309  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
310  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
312  *smolyak_basis, quad);
313 
315  smolyak_basis);
316  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
317 
318  // Test grid
320  sm_op, quad_op, setup.rtol, setup.atol, out);
321 
322  // Test apply
323  success = success &&
325  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
326 
327  // Test transpose apply
328  success = success &&
330  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
331 
332  // Test discrete orthogonality
333  success = success &&
335  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
336  }
337 
338  TEUCHOS_UNIT_TEST( Smolyak, GaussPatterson ) {
339  // Build isotropic Smolyak basis of dimension d and order p with
340  // exponential Gauss-Patterson growth
342  for (ordinal_type i=0; i<setup.d; i++)
346 
347  // Build corresponding pseudospectral operator
349  *smolyak_basis, true);
350 
351  // Generate sparse grids using original approach
352  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
353  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
355  *smolyak_basis, quad);
356 
358  smolyak_basis);
359  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
360 
361  // Test grid
363  sm_op, quad_op, setup.rtol, setup.atol, out);
364 
365  // Test apply
366  success = success &&
368  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
369 
370  // Test transpose apply
371  success = success &&
373  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
374 
375  // Test discrete orthogonality
376  success = success &&
378  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
379  }
380 
381 #endif
382 
383 }
384 
385 int main( int argc, char* argv[] ) {
386  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
388 }
size_type size() const
TEUCHOS_UNIT_TEST(Coefficients, IsotropicLinear)
Legendre polynomial basis using Gauss-Patterson quadrature points.
static int runUnitTestsFromMain(int argc, char *argv[])
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
value_type quad_func2(const Teuchos::ArrayView< const value_type > &x)
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
Legendre polynomial basis.
int main(int argc, char **argv)
expr val()
An isotropic total order index set.
UnitTestSetup< ordinal_type, value_type > setup_type
bool testPseudoSpectralDiscreteOrthogonality(const basis_type &basis, const operator_type &op, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
Legendre polynomial basis using Clenshaw-Curtis quadrature points.
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
An operator for building pseudo-spectral coefficients using an arbitrary quadrature rule...
bool testPseudoSpectralPoints(const operator_type1 &op1, const operator_type2 &op2, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
bool testPseudoSpectralApply(const operator_type1 &op1, const operator_type2 &op2, const func_type1 &func1, const func_type2 &func2, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...
bool testPseudoSpectralApplyTrans(const operator_type1 &op1, const operator_type2 &op2, const func_type1 &func1, const func_type2 &func2, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
value_type quad_func1(const Teuchos::ArrayView< const value_type > &x)