Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_HouseTriDiagUnitTest.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 
14 
15 #include "Stokhos.hpp"
18 
19 // Quadrature functor to be passed into quadrature expansion for mapping
20 // from Lanczos basis back to original PCE
24  pce(pce_), basis(basis_), vec(1) {}
25 
26  double operator() (const double& a) const {
27  vec[0] = a;
28  return pce.evaluate(vec);
29  }
33 };
34 
35 // Class encapsulating setup of the Lanczos basis for a given PCE
36 // u = Func(x), where Func is specified by a template-parameter
37 template <typename Func>
39  typedef typename Func::OrdinalType OrdinalType;
40  typedef typename Func::ValueType ValueType;
42  Func func;
51 
52  Lanczos_PCE_Setup(bool normalize) :
53  func()
54  {
55  rtol = 1e-8;
56  atol = 1e-12;
57  const OrdinalType d = 3;
58  const OrdinalType p = 5;
59 
60  // Create product basis
62  for (OrdinalType i=0; i<d; i++)
63  bases[i] =
65  basis =
67 
68  // Create approximation
69  sz = basis->size();
71  for (OrdinalType i=0; i<d; i++)
72  x.term(i, 1) = 1.0;
73 
74  // Tensor product quadrature
77 
78  // Triple product tensor
81 
82  // Quadrature expansion
84 
85  // Compute PCE via quadrature expansion
86  u.reset(basis);
87  v.reset(basis);
88  func.eval(*exp, x, u);
89  exp->times(v,u,u);
90 
91  // Compute Lanczos basis
92  st_bases.resize(1);
95  p, u, *Cijk));
97 
98  st_basis =
100  st_sz = st_basis->size();
105 
106  // Tensor product quadrature
107  st_quad =
109 
110  // Triple product tensor
113 
114  // Quadrature expansion
116  st_Cijk,
117  st_quad);
118 
119  st_exp.times(v_st, u_st, u_st);
120  }
121 
122 };
123 
124 #define LANCZOS_UNIT_TESTS(BASENAME, TAG, FUNC, NORMALIZE) \
125 namespace BASENAME ## TAG { \
126  \
127  Lanczos_PCE_Setup< FUNC<int,double> > setup(NORMALIZE); \
128  \
129  TEUCHOS_UNIT_TEST( BASENAME, TAG ## Map ) { \
130  Stokhos::OrthogPolyApprox<int,double> u2(setup.basis); \
131  setup.st_1d_proj_basis->transformCoeffsFromHouse( \
132  setup.u_st.coeff(), \
133  u2.coeff()); \
134  success = Stokhos::comparePCEs(setup.u, "u", u2, "u2", \
135  setup.rtol, setup.atol, out); \
136  } \
137  \
138  TEUCHOS_UNIT_TEST( BASENAME, TAG ## Orthog ) { \
139  const Teuchos::Array<double>& norms = \
140  setup.st_bases[0]->norm_squared(); \
141  const Teuchos::Array<double>& weights = \
142  setup.st_quad->getQuadWeights(); \
143  const Teuchos::Array< Teuchos::Array<double> >& values = \
144  setup.st_quad->getBasisAtQuadPoints(); \
145  Teuchos::SerialDenseMatrix<int,double> mat(setup.st_sz, \
146  setup.st_sz); \
147  for (int i=0; i<setup.st_sz; i++) { \
148  for (int j=0; j<setup.st_sz; j++) { \
149  for (unsigned int k=0; k<weights.size(); k++) \
150  mat(i,j) += weights[k]*values[k][i]*values[k][j]; \
151  mat(i,j) /= std::sqrt(norms[i]*norms[j]); \
152  } \
153  mat(i,i) -= 1.0; \
154  } \
155  success = mat.normInf() < setup.atol; \
156  if (!success) { \
157  out << "\n Error, mat.normInf() < atol = " << mat.normInf() \
158  << " < " << setup.atol << ": failed!\n"; \
159  out << "mat = " << mat << std::endl; \
160  } \
161  } \
162  \
163  TEUCHOS_UNIT_TEST( BASENAME, TAG ## PCE ) { \
164  Stokhos::OrthogPolyApprox<int,double> v2(setup.basis); \
165  lanczos_pce_quad_func quad_func(setup.v_st, *setup.st_basis); \
166  setup.exp->unary_op(quad_func, v2, setup.u); \
167  success = comparePCEs(setup.v, "v", v2, "v2", setup.rtol, \
168  setup.atol, out); \
169  } \
170  \
171  TEUCHOS_UNIT_TEST( BASENAME, TAG ## Mean ) { \
172  success = Teuchos::testRelErr( \
173  "v.mean()", setup.v.mean(), \
174  "v_st.mean()", setup.v_st.mean(), \
175  "rtol", setup.rtol, \
176  "rtol", setup.rtol, \
177  Teuchos::Ptr<std::ostream>(out.getOStream().get())); \
178  \
179  } \
180  \
181  TEUCHOS_UNIT_TEST( BASENAME, TAG ## StandardDeviation ) { \
182  success = Teuchos::testRelErr( \
183  "v.standard_deviation()", \
184  setup.v.standard_deviation(), \
185  "v_st.standard_devaition()", \
186  setup.v_st.standard_deviation(), \
187  "rtol", 1e-1, \
188  "rtol", 1e-1, \
189  Teuchos::Ptr<std::ostream>(out.getOStream().get())); \
190  } \
191  \
192 }
193 
194 //
195 // Lanczos tests based on expansion of u = cos(x) where x is a U([-1,1])
196 // random variable
197 //
198 template <typename Ordinal_Type, typename Value_Type>
200  typedef Ordinal_Type OrdinalType;
201  typedef Value_Type ValueType;
202  static const bool is_even = true;
203  void
207  exp.cos(u,x);
208  }
209 };
210 LANCZOS_UNIT_TESTS(Stokhos_LanczosPCEBasis_Proj, Cos, Lanczos_Cos_Func,
211  false)
212 
213 //
214 // Lanczos tests based on expansion of u = sin(x) where x is a U([-1,1])
215 // random variable
216 //
217 template <typename Ordinal_Type, typename Value_Type>
219  typedef Ordinal_Type OrdinalType;
220  typedef Value_Type ValueType;
221  static const bool is_even = false;
222  void
226  exp.sin(u,x);
227  }
228 };
229 LANCZOS_UNIT_TESTS(Stokhos_LanczosPCEBasis_Proj, Sin, Lanczos_Sin_Func,
230  false)
231 
232 //
233 // Lanczos tests based on expansion of u = exp(x) where x is a U([-1,1])
234 // random variable. For this test we don't use the PCE quad points and
235 // instead use those generated for the Lanczos basis
236 //
237 template <typename Ordinal_Type, typename Value_Type>
239  typedef Ordinal_Type OrdinalType;
240  typedef Value_Type ValueType;
241  static const bool is_even = false;
242  void
246  exp.exp(u,x);
247  }
248 };
249 LANCZOS_UNIT_TESTS(Stokhos_LanczosPCEBasis_Proj, Exp, Lanczos_Exp_Func,
250  false)
251 
252 int main( int argc, char* argv[] ) {
253  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
255 }
void sin(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
const Stokhos::OrthogPolyApprox< int, double > & pce
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)
lanczos_pce_quad_func(const Stokhos::OrthogPolyApprox< int, double > &pce_, const Stokhos::OrthogPolyBasis< int, double > &basis_)
value_type evaluate(const Teuchos::Array< value_type > &point) const
Evaluate polynomial approximation at a point.
Teuchos::RCP< const Stokhos::HouseTriDiagPCEBasis< OrdinalType, ValueType > > st_1d_proj_basis
void exp(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
void cos(OrthogPolyApprox< ordinal_type, value_type, node_type > &c, const OrthogPolyApprox< ordinal_type, value_type, node_type > &a)
value_type getNewCoeffs(ordinal_type i) const
Get new coefficients in this new basis.
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > v_st
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > u
void eval(Stokhos::QuadOrthogPolyExpansion< OrdinalType, ValueType > &exp, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &x, Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &u)
void eval(Stokhos::QuadOrthogPolyExpansion< OrdinalType, ValueType > &exp, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &x, Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &u)
Teuchos::RCP< const Stokhos::CompletePolynomialBasis< OrdinalType, ValueType > > st_basis
static int runUnitTestsFromMain(int argc, char *argv[])
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void reset(const Teuchos::RCP< const Stokhos::OrthogPolyBasis< ordinal_type, value_type > > &new_basis, ordinal_type sz=0)
Reset to a new basis.
Teuchos::RCP< const Stokhos::Quadrature< OrdinalType, ValueType > > st_quad
const Stokhos::OrthogPolyBasis< int, double > & basis
virtual Teuchos::RCP< Stokhos::Sparse3Tensor< ordinal_type, value_type > > computeTripleProductTensor() const
Compute triple product tensor.
void resize(size_type new_size, const value_type &x=value_type())
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > u_st
KOKKOS_INLINE_FUNCTION PCE< Storage > exp(const PCE< Storage > &a)
int main(int argc, char **argv)
Teuchos::RCP< Stokhos::QuadOrthogPolyExpansion< OrdinalType, ValueType > > exp
virtual ordinal_type size() const
Return total size of basis.
Teuchos::RCP< const Stokhos::CompletePolynomialBasis< OrdinalType, ValueType > > basis
void eval(Stokhos::QuadOrthogPolyExpansion< OrdinalType, ValueType > &exp, const Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &x, Stokhos::OrthogPolyApprox< OrdinalType, ValueType > &u)
Teuchos::Array< Teuchos::RCP< const Stokhos::OneDOrthogPolyBasis< OrdinalType, ValueType > > > st_bases
double operator()(const double &a) const
#define LANCZOS_UNIT_TESTS(BASENAME, TAG, FUNC, NORMALIZE)
Defines quadrature for a tensor product basis by tensor products of 1-D quadrature rules...
Stokhos::OrthogPolyApprox< OrdinalType, ValueType > v
reference term(ordinal_type dimension, ordinal_type order)
Get coefficient term for given dimension and order.