Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
stieltjes_example5.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 #include "Stokhos_Sacado.hpp"
15 
17 #include "Stokhos_SDMUtils.hpp"
18 
20 
21 template <typename ordinal_type, typename scalar_type>
22 void
25  double shift,
27 {
28  ordinal_type n = A.numCols();
30  s.putScalar(shift);
31  y.assign(x);
32  y -= s;
33  f.multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, 1.0, A, y, 0.0);
34 }
35 
36 template <typename ordinal_type, typename scalar_type>
39 {
40  return std::exp(-f.dot(f));
41  //return -f.dot(f);
42 }
43 
50 
51 // measure transformation approaches
53 const int num_mt_method = 3;
56 const char *mt_method_names[] = { "Stieltjes", "Lanczos", "Gram-Schmidt" };
57 
58 double rel_err(double a, double b) {
59  return std::abs(a-b)/std::abs(b);
60 }
61 
62 int main(int argc, char **argv)
63 {
64  try {
65 
66  // Setup command line options
68  CLP.setDocString(
69  "This example runs a Stieltjes-based dimension reduction example.\n");
70 
71  int n = 4;
72  CLP.setOption("n", &n, "Number of random variables");
73 
74  int m = 2;
75  CLP.setOption("m", &m, "Number of intermediate variables");
76 
77  int pmin = 1;
78  CLP.setOption("pmin", &pmin, "Starting expansion order");
79 
80  int pmax = 7;
81  CLP.setOption("pmax", &pmax, "Final expansion order");
82 
83  MT_METHOD mt_method = MT_LANCZOS;
84  CLP.setOption("mt_method", &mt_method,
86  "Measure transformation method");
87 
88  bool normalize = false;
89  CLP.setOption("normalize", "unnormalize", &normalize,
90  "Normalize PC basis");
91 
92  bool project_integrals = false;
93  CLP.setOption("project", "no_project", &project_integrals,
94  "Project integrals");
95 
96 #ifdef HAVE_STOKHOS_DAKOTA
97  bool sparse_grid = true;
98  CLP.setOption("sparse", "tensor", &sparse_grid,
99  "Use Sparse or Tensor grid");
100 #else
101  bool sparse_grid = false;
102 #endif
103 
104  double rank_threshold = 1.0e-12;
105  CLP.setOption("rank_threshold", &rank_threshold, "Rank threshold");
106 
107  double reduction_tolerance = 1.0e-12;
108  CLP.setOption("reduction_tolerance", &reduction_tolerance, "Quadrature reduction tolerance");
109 
110  bool verbose = false;
111  CLP.setOption("verbose", "quient", &verbose,
112  "Verbose output");
113 
114  CLP.parse( argc, argv );
115 
116  std::cout << "Summary of command line options:" << std::endl
117  << "\tm = " << m << std::endl
118  << "\tn = " << n << std::endl
119  << "\tmt_method = " << mt_method_names[mt_method] << std::endl
120  << "\tnormalize = " << normalize << std::endl
121  << "\tproject = " << project_integrals << std::endl
122  << "\tsparse = " << sparse_grid << std::endl
123  << "\trank_threshold = " << rank_threshold << std::endl
124  << "\treduction_tolerance = " << reduction_tolerance << std::endl;
125 
126  int np = pmax-pmin+1;
127 
128  Teuchos::Array<double> mean(np), mean_st(np), std_dev(np), std_dev_st(np);
129  Teuchos::Array<double> pt(np), pt_st(np);
130 
132 
133  Teuchos::Array<double> eval_pt(n, 0.56789);
134  double pt_true = 0.0;
135 
136  SDM0 B0(n,n), Q0, R0;
137  B0.random();
138  Teuchos::Array<double> w(n, 1.0);
139  Stokhos::QR_MGS(n, B0, w, Q0, R0);
140  SDM0 A0(Teuchos::View, B0, m, n);
141 
142  // Loop over orders
143  int idx = 0;
144  for (int p=pmin; p<=pmax; p++) {
145 
146  std::cout << "p = " << p;
147 
148  // Create product basis
149  for (int i=0; i<n; i++)
150  bases[i] = Teuchos::rcp(new basis_type(p, normalize));
153  std::cout << ", basis sz = " << basis->size();
154 
155  // Quadrature
157 #ifdef HAVE_STOKHOS_DAKOTA
158  if (sparse_grid)
159  quad =
160  Teuchos::rcp(new Stokhos::SparseGridQuadrature<int,double>(basis, p));
161 #endif
162  if (!sparse_grid)
163  quad =
165  std::cout << ", quad sz = " << quad->size();
166 
167  // Triple product tensor
169  basis->computeTripleProductTensor();
170 
171  // Quadrature expansion
174 
175  // Create approximation
176  SDV x(n);
177  for (int i=0; i<n; i++) {
178  x[i].copyForWrite();
179  x[i].reset(quad_exp);
180  if (normalize)
181  x[i].term(i, 1) = 1.0/std::sqrt(3);
182  else
183  x[i].term(i, 1) = 1.0;
184  }
185 
186  // Create matrix
187  SDM A(m,n);
188  for (int i=0; i<m; i++)
189  for (int j=0; j<n; j++)
190  A(i,j) = A0(i,j);
191 
192  // Compute f = A*(x-1)
193  SDV f(m);
194  f_func(A, x, 1.0, f);
195 
196  //std::cout << "f = " << f << std::endl;
197 
198  // Compute g = exp(-f^T*f)
199  pce_type g = g_func(f);
200 
201  // compute true point
202  SDV0 x0(n), f0(m);
203  for (int i=0; i<n; i++)
204  x0[i] = x[i].evaluate(eval_pt);
205  f_func(A0, x0, 1.0, f0);
206  pt_true = g_func(f0);
207 
208  // Compute reduced basis
209  Teuchos::ParameterList params;
210  params.set("Verbose", verbose);
211  if (mt_method == MT_GRAM_SCHMIDT)
212  params.set("Reduced Basis Method", "Monomial Proj Gram-Schmidt");
213  else if (mt_method == MT_LANCZOS)
214  params.set("Reduced Basis Method", "Product Lanczos");
215  else if (mt_method == MT_STIELTJES) {
216  params.set("Reduced Basis Method", "Product Lanczos");
217  params.set("Use Old Stieltjes Method", true);
218  }
219  params.set("Project", project_integrals);
220  params.set("Normalize", normalize);
221  params.set("Rank Threshold", rank_threshold);
222  Teuchos::ParameterList& red_quad_params =
223  params.sublist("Reduced Quadrature");
224  red_quad_params.set("Reduction Tolerance", reduction_tolerance);
225  red_quad_params.set("Verbose", verbose);
227  for (int i=0; i<m; i++)
228  pces[i] = f[i].getOrthogPolyApprox();
231  factory.createReducedBasis(p, pces, quad, Cijk);
233  gs_basis->getReducedQuadrature();
238  gs_exp_params->set("Use Quadrature for Times", true);
241  gs_basis,
242  gs_Cijk,
243  gs_quad,
244  gs_exp_params));
245  std::cout << ", red. basis sz = " <<gs_basis->size();
246  std::cout << ", red. quad sz = " << gs_quad->size();
247 
248  SDV f_st(m);
249  for (int i=0; i<m; i++) {
250  f_st(i).copyForWrite();
251  f_st(i).reset(st_quad_exp);
252  gs_basis->transformFromOriginalBasis(f(i).coeff(), f_st(i).coeff());
253  }
254 
255  // Compute g_st = exp(-f_st^T*f_st)
256  pce_type g_st = g_func(f_st);
257 
258  // Project g_st back to original basis
259  pce_type g2(quad_exp);
260  gs_basis->transformToOriginalBasis(g_st.coeff(), g2.coeff());
261 
262  // std::cout.precision(12);
263  // std::cout << g;
264  // std::cout << g2;
265  // std::cout << g_st;
266  mean[idx] = g.mean();
267  mean_st[idx] = g2.mean();
268  std_dev[idx] = g.standard_deviation();
269  std_dev_st[idx] = g2.standard_deviation();
270  pt[idx] = g.evaluate(eval_pt);
271  pt_st[idx] = g2.evaluate(eval_pt);
272  idx++;
273 
274  std::cout << std::endl;
275  }
276 
277  idx = 0;
278  int wi=10;
279  std::cout << "Statistical error:" << std::endl;
280  std::cout << "p "
281  << std::setw(wi) << "mean" << " "
282  << std::setw(wi) << "mean_st" << " "
283  << std::setw(wi) << "std_dev" << " "
284  << std::setw(wi) << "std_dev_st" << " "
285  << std::setw(wi) << "point" << " "
286  << std::setw(wi) << "point_st" << std::endl;
287  for (int p=pmin; p<pmax; p++) {
288  std::cout.precision(3);
289  std::cout.setf(std::ios::scientific);
290  std::cout << p << " "
291  << std::setw(wi) << rel_err(mean[idx], mean[np-1]) << " "
292  << std::setw(wi) << rel_err(mean_st[idx], mean[np-1]) << " "
293  << std::setw(wi) << rel_err(std_dev[idx], std_dev[np-1]) << " "
294  << std::setw(wi) << rel_err(std_dev_st[idx], std_dev[np-1])
295  << " "
296  << std::setw(wi) << rel_err(pt[idx], pt_true) << " "
297  << std::setw(wi) << rel_err(pt_st[idx], pt_true)
298  << std::endl;
299  idx++;
300  }
301 
302  }
303  catch (std::exception& e) {
304  std::cout << e.what() << std::endl;
305  }
306 }
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
const int num_mt_method
Sacado::ETPCE::OrthogPoly< double, Stokhos::StandardStorage< int, double > > pce_type
void QR_MGS(ordinal_type k, const Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &A, const Teuchos::Array< scalar_type > &w, Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &Q, Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &R)
Compute thin QR using modified Gram-Schmidt.
Teuchos::SerialDenseVector< int, double > SDV0
Generate a basis from a given set of PCE expansions that is orthogonal with respect to the product me...
ScalarType dot(const SerialDenseVector< OrdinalType, ScalarType > &x) const
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
RCP< ParameterList > sublist(const RCP< ParameterList > &paramList, const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Teuchos::SerialDenseMatrix< int, pce_type > SDM
Stokhos::LegendreBasis< int, double > basis_type
Teuchos::SerialDenseVector< int, pce_type > SDV
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
EParseCommandLineReturn parse(int argc, char *argv[], std::ostream *errout=&std::cerr) const
void f_func(const Teuchos::SerialDenseMatrix< ordinal_type, scalar_type > &A, const Teuchos::SerialDenseVector< ordinal_type, scalar_type > &x, double shift, Teuchos::SerialDenseVector< ordinal_type, scalar_type > &f)
KOKKOS_INLINE_FUNCTION PCE< Storage > abs(const PCE< Storage > &a)
const MT_METHOD mt_method_values[]
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
virtual Teuchos::RCP< Stokhos::ReducedPCEBasis< ordinal_type, value_type > > createReducedBasis(ordinal_type p, const Teuchos::Array< Stokhos::OrthogPolyApprox< ordinal_type, value_type > > &pce, const Teuchos::RCP< const Stokhos::Quadrature< ordinal_type, value_type > > &quad, const Teuchos::RCP< const Stokhos::Sparse3Tensor< ordinal_type, value_type > > &Cijk) const
Get reduced quadrature object.
Teuchos::SerialDenseMatrix< int, double > SDM0
Legendre polynomial basis.
OrdinalType numCols() const
int main(int argc, char **argv)
ScalarType f(const Teuchos::Array< ScalarType > &x, double a, double b)
void setDocString(const char doc_string[])
double rel_err(double a, double b)
scalar_type g_func(const Teuchos::SerialDenseVector< ordinal_type, scalar_type > &f)
const char * mt_method_names[]
int n
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...
ScalarType g(const Teuchos::Array< ScalarType > &x, const ScalarType &y)