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 //
4 // Stokhos Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Eric T. Phipps (etphipp@sandia.gov).
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 // ensemble_example
43 //
44 // usage:
45 // ensemble_example
46 //
47 // output:
48 // prints the Hermite Polynomial Chaos Expansion of the simple function
49 //
50 // v = 1/(log(u)^2+1)
51 //
52 // where u = 1 + 0.1*x_1 + 0.05*x_2 + 0.01*x_3 x1,x2,x3 are zero-mean
53 // unit-variance Gaussian random variables using a pseudospectral approach
54 // for computing v via embedded ensemble propagation.
55 
56 #include "Stokhos_Sacado.hpp"
58 
59 // The function to compute the polynomial chaos expansion of,
60 // written as a template function
61 template <class ScalarType>
62 ScalarType simple_function(const ScalarType& u) {
63  ScalarType z = std::log(u);
64  return 1.0/(z*z + 1.0);
65 }
66 
67 int main(int argc, char **argv)
68 {
69  // Typename of Polynomial Chaos scalar type
70  typedef Stokhos::StandardStorage<int,double> pce_storage_type;
72 
73  // Typename of ensemble scalar type
74  const int EnsembleSize = 8;
76  typedef Sacado::MP::Vector<ensemble_storage_type> ensemble_type;
77 
78  // Short-hand for several classes used below
79  using Teuchos::Array;
80  using Teuchos::RCP;
81  using Teuchos::rcp;
85  using Stokhos::Quadrature;
89 
90  try {
91 
92  // Basis of dimension 3, order 4
93  const int d = 3;
94  const int p = 4;
95  Array< RCP<const OneDOrthogPolyBasis<int,double> > > bases(d);
96  for (int i=0; i<d; i++) {
97  bases[i] = rcp(new HermiteBasis<int,double>(p));
98  }
99  RCP<const CompletePolynomialBasis<int,double> > basis =
100  rcp(new CompletePolynomialBasis<int,double>(bases));
101 
102  // Quadrature method
103  RCP<const Quadrature<int,double> > quad =
104  rcp(new TensorProductQuadrature<int,double>(basis));
105 
106  // Triple product tensor
107  RCP<Sparse3Tensor<int,double> > Cijk =
108  basis->computeTripleProductTensor();
109 
110  // Expansion method
111  RCP<QuadOrthogPolyExpansion<int,double> > expn =
112  rcp(new QuadOrthogPolyExpansion<int,double>(basis, Cijk, quad));
113 
114  // Polynomial expansion of u
115  pce_type u(expn);
116  u.term(0,0) = 1.0; // zeroth order term
117  u.term(0,1) = 0.1; // first order term for dimension 0
118  u.term(1,1) = 0.05; // first order term for dimension 1
119  u.term(2,1) = 0.01; // first order term for dimension 2
120 
121  //
122  // Compute PCE expansion of function using NISP with ensemble propagation
123  //
124 
125  // Extract quadrature data
126  const int pce_size = basis->size();
127  const int num_quad_points = quad->size();
128  const Array<double>& quad_weights = quad->getQuadWeights();
129  const Array< Array<double> >& quad_points = quad->getQuadPoints();
130  const Array< Array<double> >& quad_values = quad->getBasisAtQuadPoints();
131 
132  // Loop over quadrature points in blocks of size EnsembleSize
133  pce_type v(expn);
134  ensemble_type u_ensemble;
135  for (int qp_block=0; qp_block<num_quad_points; qp_block+=EnsembleSize) {
136  const int qp_sz = qp_block+EnsembleSize <= num_quad_points ?
137  EnsembleSize : num_quad_points-qp_block;
138 
139  // Evaluate u at each quadrature point
140  for (int qp=0; qp<qp_sz; ++qp)
141  u_ensemble.fastAccessCoeff(qp) =
142  u.evaluate(quad_points[qp_block+qp], quad_values[qp_block+qp]);
143  for (int qp=qp_sz; qp<EnsembleSize; ++qp)
144  u_ensemble.fastAccessCoeff(qp) = u_ensemble.fastAccessCoeff(qp_sz-1);
145 
146  // Evaluate function at each quadrature point
147  ensemble_type v_ensemble = simple_function(u_ensemble);
148 
149  // Sum results into PCE integral
150  for (int pc=0; pc<pce_size; ++pc) {
151  const double inv_nrm_sq = 1.0 / basis->norm_squared(pc);
152  for (int qp=0; qp<qp_sz; ++qp) {
153  const double w = quad_weights[qp_block+qp];
154  const double psi = quad_values[qp_block+qp][pc];
155  v.fastAccessCoeff(pc) +=
156  inv_nrm_sq * w * v_ensemble.fastAccessCoeff(qp) * psi;
157  }
158  }
159  }
160 
161  // Print u and v
162  std::cout << "\tu = ";
163  u.print(std::cout);
164  std::cout << "\tv = ";
165  v.print(std::cout);
166 
167  // Compute moments
168  double mean = v.mean();
169  double std_dev = v.standard_deviation();
170 
171  // Evaluate PCE and function at a point = 0.25 in each dimension
173  for (int i=0; i<d; i++)
174  pt[i] = 0.25;
175  double up = u.evaluate(pt);
176  double vp = simple_function(up);
177  double vp2 = v.evaluate(pt);
178 
179  // Print results
180  std::cout << "\tv mean = " << mean << std::endl;
181  std::cout << "\tv std. dev. = " << std_dev << std::endl;
182  std::cout << "\tv(0.25) (true) = " << vp << std::endl;
183  std::cout << "\tv(0.25) (pce) = " << vp2 << std::endl;
184  }
185  catch (std::exception& e) {
186  std::cout << e.what() << std::endl;
187  }
188 }
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...