Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gram_schmidt_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"
17 
19 
20 struct pce_quad_func {
24  pce(pce_), basis(basis_), vec(2) {}
25 
26  double operator() (const double& a, const double& b) const {
27  vec[0] = a;
28  vec[1] = b;
29  return pce.evaluate(vec);
30  }
34 };
35 
36 int main(int argc, char **argv)
37 {
38  try {
39 
40  const unsigned int d = 2;
41  const unsigned int p = 5;
42 
43  // Create product basis
45  for (unsigned int i=0; i<d; i++)
46  bases[i] = Teuchos::rcp(new basis_type(p));
49 
50  // Create approximation
51  Stokhos::OrthogPolyApprox<int,double> x(basis), u(basis), v(basis),
52  w(basis), w2(basis), w3(basis);
53  for (unsigned int i=0; i<d; i++) {
54  x.term(i, 1) = 1.0;
55  }
56 
57  // Tensor product quadrature
60 
61  // Triple product tensor
63  basis->computeTripleProductTensor();
64 
65  // Quadrature expansion
66  Stokhos::QuadOrthogPolyExpansion<int,double> quad_exp(basis, Cijk, quad);
67 
68  // Compute PCE via quadrature expansion
69  quad_exp.sin(u,x);
70  quad_exp.exp(v,x);
71  quad_exp.times(w,v,u);
72 
73  // Compute tensor product Stieltjes basis for u and v
75  st_bases[0] =
77  p, Teuchos::rcp(&u,false), quad, true));
78  st_bases[1] =
80  p, Teuchos::rcp(&v,false), quad, true));
83  Stokhos::OrthogPolyApprox<int,double> u_st(st_basis), v_st(st_basis),
84  w_st(st_basis);
85  u_st.term(0, 0) = u.mean();
86  u_st.term(0, 1) = 1.0;
87  v_st.term(0, 0) = v.mean();
88  v_st.term(1, 1) = 1.0;
89 
90  // Use Gram-Schmidt to orthogonalize Stieltjes basis of u and v
91  Teuchos::Array<double> st_points_0;
92  Teuchos::Array<double> st_weights_0;
94  st_bases[0]->getQuadPoints(p+1, st_points_0, st_weights_0, st_values_0);
95  Teuchos::Array<double> st_points_1;
96  Teuchos::Array<double> st_weights_1;
98  st_bases[1]->getQuadPoints(p+1, st_points_1, st_weights_1, st_values_1);
99  Teuchos::Array< Teuchos::Array<double> > st_points(st_points_0.size());
100  for (int i=0; i<st_points_0.size(); i++) {
101  st_points[i].resize(2);
102  st_points[i][0] = st_points_0[i];
103  st_points[i][1] = st_points_1[i];
104  }
105  Teuchos::Array<double> st_weights = st_weights_0;
106 
109  st_points,
110  st_weights,
111  1e-15));
112 
113  // Create quadrature for Gram-Schmidt basis using quad points and
114  // and weights from original basis mapped to Stieljtes basis
116  Teuchos::rcp(&st_points,false);
118  Teuchos::rcp(&st_weights,false);
121  points,
122  weights));
123 
124  // Triple product tensor
126  gs_basis->computeTripleProductTensor();
127 
128  // Gram-Schmidt quadrature expansion
129  Stokhos::QuadOrthogPolyExpansion<int,double> gs_quad_exp(gs_basis,
130  gs_Cijk,
131  gs_quad);
132 
133  Stokhos::OrthogPolyApprox<int,double> u_gs(gs_basis), v_gs(gs_basis),
134  w_gs(gs_basis);
135 
136  // Map expansion in Stieltjes basis to Gram-Schmidt basis
137  gs_basis->transformCoeffs(u_st.coeff(), u_gs.coeff());
138  gs_basis->transformCoeffs(v_st.coeff(), v_gs.coeff());
139 
140  // Compute w_gs = u_gs*v_gs in Gram-Schmidt basis
141  gs_quad_exp.times(w_gs, u_gs, v_gs);
142 
143  // Project w_gs back to original basis
144  pce_quad_func gs_func(w_gs, *gs_basis);
145  quad_exp.binary_op(gs_func, w2, u, v);
146 
147  // Triple product tensor
149  st_basis->computeTripleProductTensor();
150 
151  // Stieltjes quadrature expansion
154  points,
155  weights));
157  st_Cijk,
158  st_quad);
159 
160  // Compute w_st = u_st*v_st in Stieltjes basis
161  st_quad_exp.times(w_st, u_st, v_st);
162 
163  // Project w_st back to original basis
164  pce_quad_func st_func(w_st, *st_basis);
165  quad_exp.binary_op(st_func, w3, u, v);
166 
167  std::cout.precision(12);
168  std::cout << "w = " << std::endl << w;
169  std::cout << "w2 = " << std::endl << w2;
170  std::cout << "w3 = " << std::endl << w3;
171  std::cout << "w_gs = " << std::endl << w_gs;
172  std::cout << "w_st = " << std::endl << w_st;
173 
174  std::cout.setf(std::ios::scientific);
175  std::cout << "w.mean() = " << w.mean() << std::endl
176  << "w2.mean() = " << w2.mean() << std::endl
177  << "w3.mean() = " << w3.mean() << std::endl
178  << "w_gs.mean() = " << w_gs.mean() << std::endl
179  << "w_st.mean() = " << w_st.mean() << std::endl
180  << "w.std_dev() = " << w.standard_deviation() << std::endl
181  << "w2.std_dev() = " << w2.standard_deviation() << std::endl
182  << "w3.std_dev() = " << w3.standard_deviation() << std::endl
183  << "w_gs.std_dev() = " << w_gs.standard_deviation() << std::endl
184  << "w_st.std_dev() = " << w_st.standard_deviation() << std::endl;
185  }
186  catch (std::exception& e) {
187  std::cout << e.what() << std::endl;
188  }
189 }
void binary_op(const FuncT &func, OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
Nonlinear binary function.
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void times(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a, const OrthogPolyApprox< ordinal_type, value_type, node_type > &b)
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)
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
Generates three-term recurrence using the Discretized Stieltjes procedure applied to a polynomial cha...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const Stokhos::OrthogPolyBasis< int, double > & basis
Transforms a non-orthogonal multivariate basis to an orthogonal one using the Gram-Schmit procedure...
double operator()(const double &a, const double &b) const
const Stokhos::OrthogPolyApprox< int, double > & pce
void resize(size_type new_size, const value_type &x=value_type())
value_type mean() const
Compute mean of expansion.
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...