Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp
Go to the documentation of this file.
1 // $Id$
2 // $Source$
3 // @HEADER
4 // ***********************************************************************
5 //
6 // Stokhos Package
7 // Copyright (2009) Sandia Corporation
8 //
9 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
10 // license for use of this work by or on behalf of the U.S. Government.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
40 //
41 // ***********************************************************************
42 // @HEADER
43 
48 #include "Teuchos_ArrayView.hpp"
49 
50 #include "Stokhos.hpp"
52 
53 namespace SmolyakBasisUtilsUnitTest {
54 
55  // Common setup for unit tests
56  template <typename ordinal_type, typename value_type>
57  struct UnitTestSetup {
61 
63  rtol = 1e-12;
64  atol = 1e-12;
65  apply_rtol = 1e-2;
66  apply_atol = 1e-3;
67  p = 4;
68  d = 2;
69  }
70 
71  };
72 
73  typedef int ordinal_type;
74  typedef double value_type;
77 
78  // Function for testing quadratures
80  value_type val = 0.0;
81  for (int i=0; i<x.size(); i++)
82  val += x[i];
83  return std::exp(val);
84  }
85 
86  // Function for testing quadratures
88  value_type val = 0.0;
89  for (int i=0; i<x.size(); i++)
90  val += x[i];
91  return std::sin(val);
92  }
93 
94  TEUCHOS_UNIT_TEST( Direct, Linear ) {
95  // Build isotropic Smolyak basis of dimension d and order p with
96  // linear growth
98  for (ordinal_type i=0; i<setup.d; i++)
102 
103  // Build corresponding pseudospectral operator
105  *smolyak_basis, false);
106 
107  // Generate sparse grids using original approach
108  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
109  smolyak_basis, setup.p, 1e-12, Pecos::SLOW_RESTRICTED_GROWTH);
111  *smolyak_basis, quad);
112 
113  // Test grid
115  sm_op, quad_op, setup.rtol, setup.atol, out);
116 
117  // Test apply
118  success = success &&
120  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
121 
122  // Test transpose apply
123  success = success &&
125  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
126 
127  // Test discrete orthogonality
128  success = success &&
130  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
131  }
132 
133  TEUCHOS_UNIT_TEST( Direct, ModerateLinear ) {
134  // Build isotropic Smolyak basis of dimension d and order p with
135  // moderate linear growth
137  for (ordinal_type i=0; i<setup.d; i++)
141 
142  // Build corresponding pseudospectral operator
144  *smolyak_basis, false);
145 
146  // Generate sparse grids using original approach
147  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
148  smolyak_basis, setup.p, 1e-12, Pecos::MODERATE_RESTRICTED_GROWTH);
150  *smolyak_basis, quad);
151 
152  // Test grid
154  sm_op, quad_op, setup.rtol, setup.atol, out);
155 
156  // Test apply
157  success = success &&
159  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
160 
161  // Test transpose apply
162  success = success &&
164  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
165 
166  // Direct apply will not satisfy discrete orthogonality
167  }
168 
169  TEUCHOS_UNIT_TEST( Smolyak, Linear ) {
170  // Build isotropic Smolyak basis of dimension d and order p with
171  // linear growth
173  for (ordinal_type i=0; i<setup.d; i++)
177 
178  // Build corresponding pseudospectral operator
180  *smolyak_basis, true);
181 
182  // Generate sparse grids using original approach
183  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
184  smolyak_basis, setup.p, 1e-12, Pecos::SLOW_RESTRICTED_GROWTH);
186  *smolyak_basis, quad);
187 
189  smolyak_basis);
190  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
191 
192  // Test grid
194  sm_op, quad_op, setup.rtol, setup.atol, out);
195 
196  // Test apply
197  success = success &&
199  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
200 
201  // Test transpose apply
202  success = success &&
204  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
205 
206  // Test discrete orthogonality
207  success = success &&
209  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
210  }
211 
212  TEUCHOS_UNIT_TEST( Smolyak, ModerateLinear ) {
213  // Build isotropic Smolyak basis of dimension d and order p with
214  // moderate linear growth
216  for (ordinal_type i=0; i<setup.d; i++)
220 
221  // Build corresponding pseudospectral operator
223  *smolyak_basis, true);
224 
225  // Generate sparse grids using original approach
226  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
227  smolyak_basis, setup.p, 1e-12, Pecos::MODERATE_RESTRICTED_GROWTH);
229  *smolyak_basis, quad);
230 
232  smolyak_basis);
233  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
234 
235  // Test grid
237  sm_op, quad_op, setup.rtol, setup.atol, out);
238 
239  // Test apply
240  success = success &&
242  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
243 
244  // Test transpose apply
245  success = success &&
247  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
248 
249  // Test discrete orthogonality
250  success = success &&
252  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
253  }
254 
255 #ifdef HAVE_STOKHOS_DAKOTA
256 
257  TEUCHOS_UNIT_TEST( Direct, ClenshawCurtis ) {
258  // Build isotropic Smolyak basis of dimension d and order p with
259  // exponential growth
261  for (ordinal_type i=0; i<setup.d; i++)
265 
266  // Build corresponding pseudospectral operator
268  *smolyak_basis, false);
269 
270  // Generate sparse grids using original approach
271  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
272  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
274  *smolyak_basis, quad);
275 
276  // Test grid
278  sm_op, quad_op, setup.rtol, setup.atol, out);
279 
280  // Test apply
281  success = success &&
283  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
284 
285  // Test transpose apply
286  success = success &&
288  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
289 
290  // Direct apply will not satisfy discrete orthogonality
291  }
292 
293  TEUCHOS_UNIT_TEST( Direct, GaussPatterson ) {
294  // Build isotropic Smolyak basis of dimension d and order p with
295  // exponential Gauss-Patterson growth
297  for (ordinal_type i=0; i<setup.d; i++)
301 
302  // Build corresponding pseudospectral operator
304  *smolyak_basis, false);
305 
306  // Generate sparse grids using original approach
307  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
308  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
310  *smolyak_basis, quad);
311 
312  // Test grid
314  sm_op, quad_op, setup.rtol, setup.atol, out);
315 
316  // Test apply
317  success = success &&
319  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
320 
321  // Test transpose apply
322  success = success &&
324  sm_op, quad_op, quad_func1, quad_func2, setup.rtol, setup.atol, out);
325 
326  // Direct apply will not satisfy discrete orthogonality
327  }
328 
329  TEUCHOS_UNIT_TEST( Smolyak, ClenshawCurtis ) {
330  // Build isotropic Smolyak basis of dimension d and order p with
331  // exponential growth
333  for (ordinal_type i=0; i<setup.d; i++)
337 
338  // Build corresponding pseudospectral operator
340  *smolyak_basis, true);
341 
342  // Generate sparse grids using original approach
343  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
344  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
346  *smolyak_basis, quad);
347 
349  smolyak_basis);
350  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
351 
352  // Test grid
354  sm_op, quad_op, setup.rtol, setup.atol, out);
355 
356  // Test apply
357  success = success &&
359  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
360 
361  // Test transpose apply
362  success = success &&
364  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
365 
366  // Test discrete orthogonality
367  success = success &&
369  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
370  }
371 
372  TEUCHOS_UNIT_TEST( Smolyak, GaussPatterson ) {
373  // Build isotropic Smolyak basis of dimension d and order p with
374  // exponential Gauss-Patterson growth
376  for (ordinal_type i=0; i<setup.d; i++)
380 
381  // Build corresponding pseudospectral operator
383  *smolyak_basis, true);
384 
385  // Generate sparse grids using original approach
386  Stokhos::SparseGridQuadrature<ordinal_type,value_type> quad(
387  smolyak_basis, setup.p, 1e-12, Pecos::UNRESTRICTED_GROWTH);
389  *smolyak_basis, quad);
390 
392  smolyak_basis);
393  Stokhos::QuadraturePseudoSpectralOperator<ordinal_type,value_type> tp_quad_op(*smolyak_basis, tp_quad);
394 
395  // Test grid
397  sm_op, quad_op, setup.rtol, setup.atol, out);
398 
399  // Test apply
400  success = success &&
402  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
403 
404  // Test transpose apply
405  success = success &&
407  sm_op, tp_quad_op, quad_func1, quad_func2, setup.apply_rtol, setup.apply_atol, out);
408 
409  // Test discrete orthogonality
410  success = success &&
412  *smolyak_basis, sm_op, setup.rtol, setup.atol, out);
413  }
414 
415 #endif
416 
417 }
418 
419 int main( int argc, char* argv[] ) {
420  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
422 }
Stokhos::testPseudoSpectralApplyTrans
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)
Definition: Stokhos_UnitTestHelpers.hpp:465
Teuchos_TestingHelpers.hpp
SmolyakBasisUtilsUnitTest::UnitTestSetup::apply_rtol
value_type apply_rtol
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:59
Stokhos::SmolyakPseudoSpectralOperator
An operator for building pseudo-spectral coefficients using a sparse Smolyak construction.
Definition: Stokhos_SmolyakPseudoSpectralOperator.hpp:61
SmolyakBasisUtilsUnitTest::ordinal_type
int ordinal_type
Definition: Stokhos_SmolyakBasisUnitTest.cpp:70
Stokhos::QuadraturePseudoSpectralOperator
An operator for building pseudo-spectral coefficients using an arbitrary quadrature rule.
Definition: Stokhos_QuadraturePseudoSpectralOperator.hpp:60
Sacado::UQ::exp
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
Definition: Sacado_UQ_PCE_Imp.hpp:827
SmolyakBasisUtilsUnitTest::UnitTestSetup::atol
value_type atol
Definition: Stokhos_SmolyakBasisUnitTest.cpp:58
Stokhos.hpp
SmolyakBasisUtilsUnitTest::TEUCHOS_UNIT_TEST
TEUCHOS_UNIT_TEST(Coefficients, IsotropicLinear)
Definition: Stokhos_SmolyakBasisUnitTest.cpp:104
SmolyakBasisUtilsUnitTest::setup
setup_type setup
Definition: Stokhos_SmolyakBasisUnitTest.cpp:73
SmolyakBasisUtilsUnitTest::value_type
double value_type
Definition: Stokhos_SmolyakBasisUnitTest.cpp:71
Teuchos::rcp
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
SmolyakBasisUtilsUnitTest::UnitTestSetup::apply_atol
value_type apply_atol
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:59
Teuchos_ArrayView.hpp
Stokhos::SmolyakBasis
Multivariate orthogonal polynomial basis generated from a Smolyak sparse grid.
Definition: Stokhos_SmolyakBasis.hpp:62
Sacado::UQ::sin
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
Definition: Sacado_UQ_PCE_Imp.hpp:965
main
int main(int argc, char *argv[])
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:419
SmolyakBasisUtilsUnitTest::UnitTestSetup::p
ordinal_type p
Definition: Stokhos_SmolyakBasisUnitTest.cpp:59
Teuchos::ArrayView
Teuchos::GlobalMPISession
Teuchos::RCP
Stokhos::GaussPattersonLegendreBasis
Legendre polynomial basis using Gauss-Patterson quadrature points.
Definition: Stokhos_GaussPattersonLegendreBasis.hpp:57
Stokhos::TotalOrderIndexSet
An isotropic total order index set.
Definition: Stokhos_ProductBasisUtils.hpp:215
Teuchos_UnitTestHarness.hpp
Teuchos::Array
SmolyakBasisUtilsUnitTest::UnitTestSetup::d
ordinal_type d
Definition: Stokhos_SmolyakBasisUnitTest.cpp:59
Teuchos_UnitTestRepository.hpp
val
expr val()
SmolyakBasisUtilsUnitTest
Definition: Stokhos_SmolyakBasisUnitTest.cpp:53
Stokhos::TensorProductQuadrature
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules.
Definition: Stokhos_TensorProductQuadrature.hpp:58
SmolyakBasisUtilsUnitTest::quad_func2
value_type quad_func2(const Teuchos::ArrayView< const value_type > &x)
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:87
Stokhos::LegendreBasis
Legendre polynomial basis.
Definition: Stokhos_LegendreBasis.hpp:67
SmolyakBasisUtilsUnitTest::UnitTestSetup::UnitTestSetup
UnitTestSetup()
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:62
UnitTestSetup
Definition: Stokhos_SacadoMPVectorCommTests.cpp:71
Teuchos_GlobalMPISession.hpp
Stokhos_UnitTestHelpers.hpp
SmolyakBasisUtilsUnitTest::UnitTestSetup::rtol
value_type rtol
Definition: Stokhos_SmolyakBasisUnitTest.cpp:58
SmolyakBasisUtilsUnitTest::setup_type
UnitTestSetup< ordinal_type, value_type > setup_type
Definition: Stokhos_SmolyakBasisUnitTest.cpp:72
Stokhos::MODERATE_GROWTH
Definition: Stokhos_RecurrenceBasis.hpp:52
cusp::detail::device::x
const IndexType const IndexType const IndexType const IndexType const ValueType const ValueType * x
Definition: csr_vector.h:260
Stokhos::testPseudoSpectralPoints
bool testPseudoSpectralPoints(const operator_type1 &op1, const operator_type2 &op2, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
Definition: Stokhos_UnitTestHelpers.hpp:362
SmolyakBasisUtilsUnitTest::quad_func1
value_type quad_func1(const Teuchos::ArrayView< const value_type > &x)
Definition: Stokhos_SmolyakPseudoSpectralOperatorUnitTest.cpp:79
Teuchos::UnitTestRepository::runUnitTestsFromMain
static int runUnitTestsFromMain(int argc, char *argv[])
Stokhos::testPseudoSpectralApply
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)
Definition: Stokhos_UnitTestHelpers.hpp:404
argv
char * argv[]
Definition: Stokhos_HouseTriDiagUnitTest.cpp:286
Stokhos::testPseudoSpectralDiscreteOrthogonality
bool testPseudoSpectralDiscreteOrthogonality(const basis_type &basis, const operator_type &op, const scalar_type &rel_tol, const scalar_type &abs_tol, Teuchos::FancyOStream &out)
Definition: Stokhos_UnitTestHelpers.hpp:524
Stokhos::ClenshawCurtisLegendreBasis
Legendre polynomial basis using Clenshaw-Curtis quadrature points.
Definition: Stokhos_ClenshawCurtisLegendreBasis.hpp:57