Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
stieltjes_example3.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 #include <iostream>
43 #include <iomanip>
44 
45 #include "Stokhos.hpp"
50 
51 //typedef Stokhos::LegendreBasis<int,double> basis_type;
53 
54 struct pce_quad_func {
58  pce(pce_), basis(basis_), vec(1) {}
59 
60  double operator() (const double& a) const {
61  vec[0] = a;
62  return pce.evaluate(vec);
63  }
67 };
68 
69 struct sin_func {
71  pce(pce_) {}
72  double operator() (const Teuchos::Array<double>& x) const {
73  return std::sin(pce.evaluate(x));
74  }
76 };
77 
78 double rel_err(double a, double b) {
79  return std::abs(a-b)/std::abs(b);
80 }
81 
82 int main(int argc, char **argv)
83 {
84  try {
85 
86  const unsigned int d = 2;
87  const unsigned int pmin = 2;
88  const unsigned int pmax = 10;
89  const unsigned int np = pmax-pmin+1;
90  bool use_pce_quad_points = false;
91  bool normalize = true;
92  bool sparse_grid = true;
93 #ifndef HAVE_STOKHOS_DAKOTA
94  sparse_grid = false;
95 #endif
96  Teuchos::Array<double> mean(np), mean_st(np), std_dev(np), std_dev_st(np);
97  Teuchos::Array<double> pt(np), pt_st(np);
98 
100  Teuchos::Array<double> eval_pt(d, 0.5);
101  double pt_true;
102 
103  // Loop over orders
104  unsigned int n = 0;
105  for (unsigned int p=pmin; p<=pmax; p++) {
106 
107  std::cout << "p = " << p << std::endl;
108 
109  // Create product basis
110  for (unsigned int i=0; i<d; i++)
111  bases[i] = Teuchos::rcp(new basis_type(p));
114 
115  // Create approximation
116  Stokhos::OrthogPolyApprox<int,double> x(basis), u(basis), v(basis),
117  w(basis), w2(basis);
118  for (unsigned int i=0; i<d; i++) {
119  x.term(i, 1) = 1.0;
120  }
121 
122  double x_pt = x.evaluate(eval_pt);
123  pt_true = std::exp(std::sin(x_pt));
124 
125  // Quadrature
127 #ifdef HAVE_STOKHOS_DAKOTA
128  if (sparse_grid)
129  quad =
130  Teuchos::rcp(new Stokhos::SparseGridQuadrature<int,double>(basis, p));
131 #endif
132  if (!sparse_grid)
133  quad =
135 
136  // Triple product tensor
138  basis->computeTripleProductTensor();
139 
140  // Quadrature expansion
141  Stokhos::QuadOrthogPolyExpansion<int,double> quad_exp(basis, Cijk, quad);
142 
143  // Compute PCE via quadrature expansion
144  quad_exp.sin(u,x);
145  //quad_exp.times(u,u,10.0);
146  quad_exp.exp(w,u);
147 
148  // Compute Stieltjes basis
152  p, Teuchos::rcp(&u,false), Cijk, normalize));
153  st_bases[0] = st_1d_basis;
154 
156  st_basis =
158  //std::cout << *st_basis << std::endl;
159 
160  Stokhos::OrthogPolyApprox<int,double> u_st(st_basis), w_st(st_basis);
161  u_st.term(0, 0) = st_1d_basis->getNewCoeffs(0);
162  u_st.term(0, 1) = st_1d_basis->getNewCoeffs(1);
163 
164  // Triple product tensor
166  st_basis->computeTripleProductTensor();
167 
168  // Tensor product quadrature
170  if (!use_pce_quad_points) {
171 #ifdef HAVE_STOKHOS_DAKOTA
172  if (sparse_grid)
173  st_quad = Teuchos::rcp(new Stokhos::SparseGridQuadrature<int,double>(st_basis, p));
174 #endif
175  if (!sparse_grid)
177  }
178  else {
179  Teuchos::Array<double> st_points_0;
180  Teuchos::Array<double> st_weights_0;
182  st_bases[0]->getQuadPoints(p+1, st_points_0, st_weights_0, st_values_0);
183  Teuchos::Array<double> st_points_1;
184  Teuchos::Array<double> st_weights_1;
186  st_bases[1]->getQuadPoints(p+1, st_points_1, st_weights_1, st_values_1);
189  for (int i=0; i<st_points_0.size(); i++) {
190  (*st_points)[i].resize(2);
191  (*st_points)[i][0] = st_points_0[i];
192  (*st_points)[i][1] = st_points_1[i];
193  }
195  Teuchos::rcp(new Teuchos::Array<double>(st_weights_0));
197  st_basis;
198  st_quad =
200  st_points,
201  st_weights));
202  }
203 
204  // Quadrature expansion
205  Stokhos::QuadOrthogPolyExpansion<int,double> st_quad_exp(st_basis,
206  st_Cijk,
207  st_quad);
208 
209  // Compute w_st = u_st*v_st in Stieltjes basis
210  st_quad_exp.exp(w_st, u_st);
211 
212  // Project w_st back to original basis
213  pce_quad_func st_func(w_st, *st_basis);
214  quad_exp.unary_op(st_func, w2, u);
215 
216  // std::cout.precision(12);
217  // std::cout << w;
218  // std::cout << w2;
219  // std::cout << w_st;
220  mean[n] = w.mean();
221  mean_st[n] = w2.mean();
222  std_dev[n] = w.standard_deviation();
223  std_dev_st[n] = w2.standard_deviation();
224  pt[n] = w.evaluate(eval_pt);
225  pt_st[n] = w2.evaluate(eval_pt);
226  n++;
227  }
228 
229  n = 0;
230  int wi=10;
231  std::cout << "Statistical error:" << std::endl;
232  std::cout << "p "
233  << std::setw(wi) << "mean" << " "
234  << std::setw(wi) << "mean_st" << " "
235  << std::setw(wi) << "std_dev" << " "
236  << std::setw(wi) << "std_dev_st" << " "
237  << std::setw(wi) << "point" << " "
238  << std::setw(wi) << "point_st" << std::endl;
239  for (unsigned int p=pmin; p<pmax; p++) {
240  std::cout.precision(3);
241  std::cout.setf(std::ios::scientific);
242  std::cout << p << " "
243  << std::setw(wi) << rel_err(mean[n], mean[np-1]) << " "
244  << std::setw(wi) << rel_err(mean_st[n], mean[np-1]) << " "
245  << std::setw(wi) << rel_err(std_dev[n], std_dev[np-1]) << " "
246  << std::setw(wi) << rel_err(std_dev_st[n], std_dev[np-1])
247  << " "
248  << std::setw(wi) << rel_err(pt[n], pt_true) << " "
249  << std::setw(wi) << rel_err(pt_st[n], pt_true)
250  << std::endl;
251  n++;
252  }
253 
254  }
255  catch (std::exception& e) {
256  std::cout << e.what() << std::endl;
257  }
258 }
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
value_type evaluate(const Teuchos::Array< value_type > &point) const
Evaluate polynomial approximation at a point.
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Stokhos::OrthogPolyApprox< int, double > & pce
sin_func(const Stokhos::OrthogPolyApprox< int, double > &pce_)
pce_quad_func(const Stokhos::OrthogPolyApprox< int, double > &pce_, const Stokhos::OrthogPolyBasis< int, double > &basis_)
Teuchos::Array< double > vec
Stokhos::LegendreBasis< int, double > basis_type
value_type standard_deviation() const
Compute standard deviation of expansion.
double operator()(const Teuchos::Array< double > &x) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const Stokhos::OrthogPolyBasis< int, double > & basis
double operator()(const double &a, const double &b) const
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
const Stokhos::OrthogPolyApprox< int, double > & pce
value_type mean() const
Compute mean of expansion.
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
int main(int argc, char **argv)
Legendre polynomial basis using Clenshaw-Curtis quadrature points.
size_type size() const
double rel_err(double a, double b)
KOKKOS_INLINE_FUNCTION PCE< Storage > sin(const PCE< Storage > &a)
int n
Generates three-term recurrence using the Lanczos procedure applied to a polynomial chaos expansion i...
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...
void unary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
Nonlinear unary function.