Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
exp_moment_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 <iostream>
11 #include <iomanip>
12 
13 #include "Stokhos.hpp"
14 
16 
17 // This example uses PC expansions for computing moments of
18 //
19 // u = exp(x1 + ... + xd)
20 //
21 // where x1, ..., xd are uniform random variables on [-1,1]. The methods
22 // are compared to computing the "true" value via very high-order quadrature.
23 // Because of the structure of the exponential, the moments can easily
24 // be computed in one dimension.
25 
26 int main(int argc, char **argv)
27 {
28  try {
29 
30  // Compute "true" 1-D mean, std. dev using quadrature
31  const unsigned int true_quad_order = 200;
32  basis_type tmp_basis(1);
33  Teuchos::Array<double> true_quad_points, true_quad_weights;
34  Teuchos::Array< Teuchos::Array<double> > true_quad_values;
35  tmp_basis.getQuadPoints(true_quad_order, true_quad_points,
36  true_quad_weights, true_quad_values);
37  double mean_1d = 0.0;
38  double sd_1d = 0.0;
39  for (unsigned int qp=0; qp<true_quad_points.size(); qp++) {
40  double t = std::exp(true_quad_points[qp]);
41  mean_1d += t*true_quad_weights[qp];
42  sd_1d += t*t*true_quad_weights[qp];
43  }
44 
45  const unsigned int dmin = 1;
46  const unsigned int dmax = 4;
47  const unsigned int pmin = 1;
48  const unsigned int pmax = 5;
49 
50  // Loop over dimensions
51  for (unsigned int d=dmin; d<=dmax; d++) {
52 
53  // compute "true" values
54  double true_mean = std::pow(mean_1d,static_cast<double>(d));
55  double true_sd = std::pow(sd_1d,static_cast<double>(d)) -
56  true_mean*true_mean;
57  true_sd = std::sqrt(true_sd);
58  std::cout.precision(12);
59  std::cout << "true mean = " << true_mean << "\t true std. dev. = "
60  << true_sd << std::endl;
61 
63 
64  // Loop over orders
65  for (unsigned int p=pmin; p<=pmax; p++) {
66 
67  // Create product basis
68  for (unsigned int i=0; i<d; i++)
69  bases[i] = Teuchos::rcp(new basis_type(p));
72 
73  // Create approximation
74  int sz = basis->size();
75  Stokhos::OrthogPolyApprox<int,double> x(basis), u(basis);
76  for (unsigned int i=0; i<d; i++) {
77  x.term(i,1) = 1.0;
78  }
79 
80  // Tensor product quadrature
83 
84  // Triple product tensor
86  basis->computeTripleProductTensor();
87 
88  // Quadrature expansion
90  quad);
91 
92  // Compute PCE via quadrature expansion
93  quad_exp.exp(u,x);
94  double mean = u.mean();
95  double sd = u.standard_deviation();
96 
97  std::cout.precision(4);
98  std::cout.setf(std::ios::scientific);
99  std::cout << "d = " << d << " p = " << p
100  << " sz = " << sz
101  << "\tmean err = "
102  << std::fabs(true_mean-mean) << "\tstd. dev. err = "
103  << std::fabs(true_sd-sd) << std::endl;
104  }
105 
106  }
107  }
108  catch (std::exception& e) {
109  std::cout << e.what() << std::endl;
110  }
111 }
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
KOKKOS_INLINE_FUNCTION PCE< Storage > fabs(const PCE< Storage > &a)
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)
virtual void getQuadPoints(ordinal_type quad_order, Teuchos::Array< value_type > &points, Teuchos::Array< value_type > &weights, Teuchos::Array< Teuchos::Array< value_type > > &values) const
Compute quadrature points, weights, and values of basis polynomials at given set of points points...
Stokhos::LegendreBasis< int, double > basis_type
value_type standard_deviation() const
Compute standard deviation of expansion.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
value_type mean() const
Compute mean of expansion.
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
Legendre polynomial basis.
int main(int argc, char **argv)
size_type size() const
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...