Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
sacado_example.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 #include "Stokhos_Sacado.hpp"
12 
13 // The function to compute the polynomial chaos expansion of,
14 // written as a template function
15 template <class ScalarType>
16 ScalarType simple_function(const ScalarType& u) {
17  ScalarType z = std::log(u);
18  return 1.0/(z*z + 1.0);
19 }
20 
21 int main(int argc, char **argv)
22 {
23  // Typename of Polynomial Chaos scalar type
26 
27  // Short-hand for several classes used below
28  using Teuchos::Array;
29  using Teuchos::RCP;
30  using Teuchos::rcp;
35  using Stokhos::Quadrature;
41 
42  try {
43 
44  // Setup command line options
46  CLP.setDocString(
47  "This example computes the PC expansion of a simple function.\n");
48  int p = 4;
49  CLP.setOption("order", &p, "Polynomial order");
50  bool sparse = false;
51  CLP.setOption("sparse", "tensor", &sparse,
52  "Use sparse grid or tensor product quadrature");
53 
54  // Parse arguments
55  CLP.parse( argc, argv );
56 
57  // Basis of dimension 3, order given by command-line option
58  const int d = 3;
59  Array< RCP<const OneDOrthogPolyBasis<int,double> > > bases(d);
60  for (int i=0; i<d; i++) {
61  bases[i] = rcp(new HermiteBasis<int,double>(p, true));
62  }
63  RCP<const CompletePolynomialBasis<int,double> > basis =
64  rcp(new CompletePolynomialBasis<int,double>(bases));
65  std::cout << "basis size = " << basis->size() << std::endl;
66 
67  // Quadrature method
68  RCP<const Quadrature<int,double> > quad;
69  if (sparse) {
70  const TotalOrderIndexSet<int> index_set(d, p);
71  quad = rcp(new SmolyakSparseGridQuadrature<int,double>(basis, index_set));
72  }
73  else {
74  quad = rcp(new TensorProductQuadrature<int,double>(basis));
75  }
76  std::cout << "quadrature size = " << quad->size() << std::endl;
77 
78  // Triple product tensor
79  RCP<Sparse3Tensor<int,double> > Cijk =
80  basis->computeTripleProductTensor();
81 
82  // Expansion method
83  RCP<QuadOrthogPolyExpansion<int,double> > expn =
84  rcp(new QuadOrthogPolyExpansion<int,double>(basis, Cijk, quad));
85 
86  // Polynomial expansion of u (note: these are coefficients in the
87  // normalized basis)
88  pce_type u(expn);
89  u.term(0,0) = 1.0; // zeroth order term
90  u.term(0,1) = 0.1; // first order term for dimension 0
91  u.term(1,1) = 0.05; // first order term for dimension 1
92  u.term(2,1) = 0.01; // first order term for dimension 2
93 
94  // Compute PCE expansion of function
95  pce_type v = simple_function(u);
96 
97  // Print u and v
98  std::cout << "\tu = ";
99  u.print(std::cout);
100  std::cout << "\tv = ";
101  v.print(std::cout);
102 
103  // Compute moments
104  double mean = v.mean();
105  double std_dev = v.standard_deviation();
106 
107  // Evaluate PCE and function at a point = 0.25 in each dimension
109  for (int i=0; i<d; i++)
110  pt[i] = 0.25;
111  double up = u.evaluate(pt);
112  double vp = simple_function(up);
113  double vp2 = v.evaluate(pt);
114 
115  // Print results
116  std::cout << "\tv mean = " << mean << std::endl;
117  std::cout << "\tv std. dev. = " << std_dev << std::endl;
118  std::cout << "\tv(0.25) (true) = " << vp << std::endl;
119  std::cout << "\tv(0.25) (pce) = " << vp2 << std::endl;
120 
121  // Check the answer
122  if (std::abs(vp - vp2) < 1e-2)
123  std::cout << "\nExample Passed!" << std::endl;
124  }
125  catch (std::exception& e) {
126  std::cout << e.what() << std::endl;
127  }
128 }
ScalarType simple_function(const ScalarType &u)
Stokhos::StandardStorage< int, double > storage_type
Hermite polynomial basis.
Sacado::ETPCE::OrthogPoly< double, Stokhos::StandardStorage< int, double > > pce_type
Data structure storing a sparse 3-tensor C(i,j,k) in a a compressed format.
Defines quadrature for a tensor product basis by Smolyak sparse grids.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract base class for quadrature methods.
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
EParseCommandLineReturn parse(int argc, char *argv[], std::ostream *errout=&std::cerr) const
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
Legendre polynomial basis.
int main(int argc, char **argv)
An isotropic total order index set.
Abstract base class for 1-D orthogonal polynomials.
void setDocString(const char doc_string[])
Orthogonal polynomial expansions based on numerical quadrature.
KOKKOS_INLINE_FUNCTION PCE< Storage > log(const PCE< Storage > &a)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...