Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_LogNormalUnitTest.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 
10 // Test the exponential of a linear Hermite expansion using an analytic
11 // closed form solution
12 
17 
18 #include "Stokhos.hpp"
20 
21 TEUCHOS_UNIT_TEST( Stokhos_LogNormal_TP, UnitTest ) {
22  // Basis of dimension 3, order 5
23  const int d = 3;
24  const int p = 5;
25  const double mean = 0.2;
26  const double std_dev = 0.1;
28  for (int i=0; i<d; i++) {
30  }
33 
34  // Quadrature method
37 
38  // Triple product tensor
40  basis->computeTripleProductTensor();
41 
42  // Expansion method
43  Stokhos::QuadOrthogPolyExpansion<int,double> expn(basis, Cijk, quad);
44 
45  // Polynomial expansions
46  Stokhos::OrthogPolyApprox<int,double> u(basis), v(basis), w(basis);
47  u.term(0,0) = mean;
48  for (int i=0; i<d; i++)
49  u.term(i,1) = std_dev / (i+1);
50 
51  // Compute expansion
52  expn.exp(v,u);
53 
54  // Compute expansion analytically
55  double w_mean = mean;
56  for (int j=0; j<d; j++)
57  w_mean += u.term(j,1)*u.term(j,1)/2.0;
58  w_mean = std::exp(w_mean);
59  for (int i=0; i<basis->size(); i++) {
60  Stokhos::MultiIndex<int> multiIndex = basis->term(i);
61  double s = 1.0;
62  for (int j=0; j<d; j++)
63  s *= std::pow(u.term(j,1), multiIndex[j]);
64  w[i] = w_mean*s/basis->norm_squared(i);
65  }
66 
67  success = Stokhos::comparePCEs(v, "quad_expansion", w, "analytic",
68  1e-12, 1e-9, out);
69 }
70 
71 #ifdef HAVE_STOKHOS_DAKOTA
72 TEUCHOS_UNIT_TEST( Stokhos_LogNormal_SG, UnitTest ) {
73  // Basis of dimension 3, order 5
74  const int d = 3;
75  const int p = 5;
76  const double mean = 0.2;
77  const double std_dev = 0.1;
79  for (int i=0; i<d; i++) {
81  }
84 
85  // Quadrature method
87  Teuchos::rcp(new Stokhos::SparseGridQuadrature<int,double>(basis, p));
88 
89  // Triple product tensor
91  basis->computeTripleProductTensor();
92 
93  // Expansion method
94  Stokhos::QuadOrthogPolyExpansion<int,double> expn(basis, Cijk, quad);
95 
96  // Polynomial expansions
97  Stokhos::OrthogPolyApprox<int,double> u(basis), v(basis), w(basis);
98  u.term(0,0) = mean;
99  for (int i=0; i<d; i++)
100  u.term(i,1) = std_dev / (i+1);
101 
102  // Compute expansion
103  expn.exp(v,u);
104 
105  // Compute expansion analytically
106  double w_mean = mean;
107  for (int j=0; j<d; j++)
108  w_mean += u.term(j,1)*u.term(j,1)/2.0;
109  w_mean = std::exp(w_mean);
110  for (int i=0; i<basis->size(); i++) {
111  Stokhos::MultiIndex<int> multiIndex = basis->term(i);
112  double s = 1.0;
113  for (int j=0; j<d; j++)
114  s *= std::pow(u.term(j,1), multiIndex[j]);
115  w[i] = w_mean*s/basis->norm_squared(i);
116  }
117 
118  success = Stokhos::comparePCEs(v, "quad_expansion", w, "analytic",
119  1e-12, 1e-9, out);
120 }
121 #endif
122 
123 int main( int argc, char* argv[] ) {
124  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
126 }
Hermite polynomial basis.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > pow(const PCE< Storage > &a, const PCE< Storage > &b)
bool comparePCEs(const PCEType &a1, const std::string &a1_name, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &a2, const std::string &a2_name, const ValueType &rel_tol, const ValueType &abs_tol, Teuchos::FancyOStream &out)
A multidimensional index.
static int runUnitTestsFromMain(int argc, char *argv[])
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
int main(int argc, char **argv)
TEUCHOS_UNIT_TEST(tAdaptivityManager, test_interface)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...