Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
ensemble_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 // ensemble_example
11 //
12 // usage:
13 // ensemble_example
14 //
15 // output:
16 // prints the Hermite Polynomial Chaos Expansion of the simple function
17 //
18 // v = 1/(log(u)^2+1)
19 //
20 // where u = 1 + 0.1*x_1 + 0.05*x_2 + 0.01*x_3 x1,x2,x3 are zero-mean
21 // unit-variance Gaussian random variables using a pseudospectral approach
22 // for computing v via embedded ensemble propagation.
23 
24 #include "Stokhos_Sacado.hpp"
26 
27 // The function to compute the polynomial chaos expansion of,
28 // written as a template function
29 template <class ScalarType>
30 ScalarType simple_function(const ScalarType& u) {
31  ScalarType z = std::log(u);
32  return 1.0/(z*z + 1.0);
33 }
34 
35 int main(int argc, char **argv)
36 {
37  // Typename of Polynomial Chaos scalar type
38  typedef Stokhos::StandardStorage<int,double> pce_storage_type;
40 
41  // Typename of ensemble scalar type
42  const int EnsembleSize = 8;
44  typedef Sacado::MP::Vector<ensemble_storage_type> ensemble_type;
45 
46  // Short-hand for several classes used below
47  using Teuchos::Array;
48  using Teuchos::RCP;
49  using Teuchos::rcp;
53  using Stokhos::Quadrature;
57 
58  try {
59 
60  // Basis of dimension 3, order 4
61  const int d = 3;
62  const int p = 4;
63  Array< RCP<const OneDOrthogPolyBasis<int,double> > > bases(d);
64  for (int i=0; i<d; i++) {
65  bases[i] = rcp(new HermiteBasis<int,double>(p));
66  }
67  RCP<const CompletePolynomialBasis<int,double> > basis =
68  rcp(new CompletePolynomialBasis<int,double>(bases));
69 
70  // Quadrature method
71  RCP<const Quadrature<int,double> > quad =
72  rcp(new TensorProductQuadrature<int,double>(basis));
73 
74  // Triple product tensor
75  RCP<Sparse3Tensor<int,double> > Cijk =
76  basis->computeTripleProductTensor();
77 
78  // Expansion method
79  RCP<QuadOrthogPolyExpansion<int,double> > expn =
80  rcp(new QuadOrthogPolyExpansion<int,double>(basis, Cijk, quad));
81 
82  // Polynomial expansion of u
83  pce_type u(expn);
84  u.term(0,0) = 1.0; // zeroth order term
85  u.term(0,1) = 0.1; // first order term for dimension 0
86  u.term(1,1) = 0.05; // first order term for dimension 1
87  u.term(2,1) = 0.01; // first order term for dimension 2
88 
89  //
90  // Compute PCE expansion of function using NISP with ensemble propagation
91  //
92 
93  // Extract quadrature data
94  const int pce_size = basis->size();
95  const int num_quad_points = quad->size();
96  const Array<double>& quad_weights = quad->getQuadWeights();
97  const Array< Array<double> >& quad_points = quad->getQuadPoints();
98  const Array< Array<double> >& quad_values = quad->getBasisAtQuadPoints();
99 
100  // Loop over quadrature points in blocks of size EnsembleSize
101  pce_type v(expn);
102  ensemble_type u_ensemble;
103  for (int qp_block=0; qp_block<num_quad_points; qp_block+=EnsembleSize) {
104  const int qp_sz = qp_block+EnsembleSize <= num_quad_points ?
105  EnsembleSize : num_quad_points-qp_block;
106 
107  // Evaluate u at each quadrature point
108  for (int qp=0; qp<qp_sz; ++qp)
109  u_ensemble.fastAccessCoeff(qp) =
110  u.evaluate(quad_points[qp_block+qp], quad_values[qp_block+qp]);
111  for (int qp=qp_sz; qp<EnsembleSize; ++qp)
112  u_ensemble.fastAccessCoeff(qp) = u_ensemble.fastAccessCoeff(qp_sz-1);
113 
114  // Evaluate function at each quadrature point
115  ensemble_type v_ensemble = simple_function(u_ensemble);
116 
117  // Sum results into PCE integral
118  for (int pc=0; pc<pce_size; ++pc) {
119  const double inv_nrm_sq = 1.0 / basis->norm_squared(pc);
120  for (int qp=0; qp<qp_sz; ++qp) {
121  const double w = quad_weights[qp_block+qp];
122  const double psi = quad_values[qp_block+qp][pc];
123  v.fastAccessCoeff(pc) +=
124  inv_nrm_sq * w * v_ensemble.fastAccessCoeff(qp) * psi;
125  }
126  }
127  }
128 
129  // Print u and v
130  std::cout << "\tu = ";
131  u.print(std::cout);
132  std::cout << "\tv = ";
133  v.print(std::cout);
134 
135  // Compute moments
136  double mean = v.mean();
137  double std_dev = v.standard_deviation();
138 
139  // Evaluate PCE and function at a point = 0.25 in each dimension
141  for (int i=0; i<d; i++)
142  pt[i] = 0.25;
143  double up = u.evaluate(pt);
144  double vp = simple_function(up);
145  double vp2 = v.evaluate(pt);
146 
147  // Print results
148  std::cout << "\tv mean = " << mean << std::endl;
149  std::cout << "\tv std. dev. = " << std_dev << std::endl;
150  std::cout << "\tv(0.25) (true) = " << vp << std::endl;
151  std::cout << "\tv(0.25) (pce) = " << vp2 << std::endl;
152  }
153  catch (std::exception& e) {
154  std::cout << e.what() << std::endl;
155  }
156 }
ScalarType simple_function(const ScalarType &u)
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.
Statically allocated storage class.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Abstract base class for quadrature methods.
Multivariate orthogonal polynomial basis generated from a total-order complete-polynomial tensor prod...
int main(int argc, char **argv)
Abstract base class for 1-D orthogonal polynomials.
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...