Stokhos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Stokhos_HermiteBasisUnitTest.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"
17 #ifdef HAVE_STOKHOS_FORUQTK
18 #include "Stokhos_gaussq.h"
19 #endif
20 
21 namespace HermiteBasisUnitTest {
22 
23  // Common setup for unit tests
24  template <typename OrdinalType, typename ValueType>
25  struct UnitTestSetup {
26  ValueType rtol, atol;
27  OrdinalType p;
29 
30  UnitTestSetup() : rtol(1e-12), atol(1e-3), p(10), basis(p) {}
31 
32  };
33 
35 
36  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Order ) {
37  int order = setup.basis.order();
38  TEUCHOS_TEST_EQUALITY(order, setup.p, out, success);
39  }
40 
41  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Size ) {
42  int size = setup.basis.size();
43  TEUCHOS_TEST_EQUALITY(size, setup.p+1, out, success);
44  }
45 
46  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Norm ) {
47  const Teuchos::Array<double>& n1 = setup.basis.norm_squared();
49  n2[0] = 1.0;
50  for (int i=1; i<=setup.p; i++)
51  n2[i] = i*n2[i-1];
52  success = Stokhos::compareArrays(n1, "n1", n2, "n2", setup.rtol, setup.atol,
53  out);
54  }
55 
56  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Norm2 ) {
59  n1[0] = setup.basis.norm_squared(0);
60  n2[0] = 1.0;
61  for (int i=1; i<=setup.p; i++) {
62  n1[i] = setup.basis.norm_squared(i);
63  n2[i] = i*n2[i-1];
64  }
65  success = Stokhos::compareArrays(n1, "n1", n2, "n2", setup.rtol, setup.atol,
66  out);
67  }
68 
69  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, QuadNorm ) {
70  const Teuchos::Array<double>& n1 = setup.basis.norm_squared();
74  setup.basis.getQuadPoints(2*setup.p, x, w, v);
75  int nqp = w.size();
76  for (int i=0; i<=setup.p; i++) {
77  n2[i] = 0;
78  for (int j=0; j<nqp; j++)
79  n2[i] += w[j]*v[j][i]*v[j][i];
80  }
81  success = Stokhos::compareArrays(n1, "n1", n2, "n2", setup.rtol, setup.atol,
82  out);
83  }
84 
85  // Tests basis is orthogonal
86  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, QuadOrthog ) {
87  const Teuchos::Array<double>& norms = setup.basis.norm_squared();
90  setup.basis.getQuadPoints(2*setup.p, x, w, v);
91  int nqp = w.size();
93  for (int i=0; i<=setup.p; i++) {
94  for (int j=0; j<=setup.p; j++) {
95  for (int k=0; k<nqp; k++)
96  mat(i,j) += w[k]*v[k][i]*v[k][j];
97  mat(i,j) /= std::sqrt(norms[i]*norms[j]);
98  }
99  mat(i,i) -= 1.0;
100  }
101  success = mat.normInf() < setup.atol;
102  if (!success) {
103  out << "\n Error, mat.normInf() < atol = " << mat.normInf()
104  << " < " << setup.atol << ": failed!\n";
105  out << "mat = " << printMat(mat) << std::endl;
106  }
107  }
108 
109  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, TripleProduct ) {
111  setup.basis.computeTripleProductTensor();
112 
115  setup.basis.getQuadPoints(3*setup.p, x, w, v);
116 
117  success = true;
118  for (int i=0; i<=setup.p; i++) {
119  for (int j=0; j<=setup.p; j++) {
120  for (int k=0; k<=setup.p; k++) {
121  double c = 0.0;
122  int nqp = w.size();
123  for (int qp=0; qp<nqp; qp++)
124  c += w[qp]*v[qp][i]*v[qp][j]*v[qp][k];
125  std::stringstream ss;
126  ss << "Cijk(" << i << "," << j << "," << k << ")";
127  success = success &&
128  Stokhos::compareValues((*Cijk)(i,j,k), ss.str(), c, "c",
129  setup.rtol, setup.atol, out);
130  }
131  }
132  }
133  }
134 
135  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, DerivDoubleProduct ) {
137  setup.basis.computeDerivDoubleProductTensor();
138 
141  setup.basis.getQuadPoints(2*setup.p, x, w, v);
142  int nqp = w.size();
143  val.resize(nqp);
144  deriv.resize(nqp);
145  for (int i=0; i<nqp; i++) {
146  val[i].resize(setup.p+1);
147  deriv[i].resize(setup.p+1);
148  setup.basis.evaluateBasesAndDerivatives(x[i], val[i], deriv[i]);
149  }
150 
151  success = true;
152  for (int i=0; i<=setup.p; i++) {
153  for (int j=0; j<=setup.p; j++) {
154  double b = 0.0;
155  for (int qp=0; qp<nqp; qp++)
156  b += w[qp]*deriv[qp][i]*val[qp][j];
157  std::stringstream ss;
158  ss << "Bij(" << i << "," << j << ")";
159  success = success &&
160  Stokhos::compareValues((*Bij)(i,j), ss.str(), b, "b",
161  setup.rtol, setup.atol, out);
162  }
163  }
164  }
165 
166  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, DerivDoubleProduct2 ) {
168  setup.basis.computeDerivDoubleProductTensor();
169  const Teuchos::Array<double>& n = setup.basis.norm_squared();
170 
171  success = true;
172  for (int i=0; i<=setup.p; i++) {
173  for (int j=0; j<=setup.p; j++) {
174  double b = 0.0;
175  if (j == i-1)
176  b = i*n[j];
177  std::stringstream ss;
178  ss << "Bij(" << i << "," << j << ")";
179  success = success &&
180  Stokhos::compareValues((*Bij)(i,j), ss.str(), b, "b",
181  setup.rtol, setup.atol, out);
182  }
183  }
184  }
185 
186  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, EvaluateBases ) {
187  double x = 0.1234;
188  Teuchos::Array<double> v1(setup.p+1), v2(setup.p+1);
189  setup.basis.evaluateBases(x, v1);
190 
191  // evaluate bases using formula
192  // He_0(x) = 1
193  // He_1(x) = x
194  // He_i(x) = x*He_{i-1}(x) - (i-1)*He_{i-2}(x), i=2,3,...
195  v2[0] = 1.0;
196  if (setup.p >= 1)
197  v2[1] = x;
198  for (int i=2; i<=setup.p; i++)
199  v2[i] = x*v2[i-1] - (i-1.0)*v2[i-2];
200  success = Stokhos::compareArrays(v1, "v1", v2, "v2", setup.rtol, setup.atol,
201  out);
202  }
203 
204  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, EvaluateBasesAndDerivatives ) {
205  double x = 0.1234;
206  Teuchos::Array<double> val1(setup.p+1), deriv1(setup.p+1),
207  val2(setup.p+1), deriv2(setup.p+1);
208  setup.basis.evaluateBasesAndDerivatives(x, val1, deriv1);
209 
210  // evaluate bases and derivatives using formula:
211  // He_0(x) = 1
212  // He_1(x) = x
213  // He_i(x) = x*He_{i-1}(x) - (i-1)*He_{i-2}(x), i=2,3,...
214  val2[0] = 1.0;
215  if (setup.p >= 1)
216  val2[1] = x;
217  for (int i=2; i<=setup.p; i++)
218  val2[i] = x*val2[i-1] - (i-1.0)*val2[i-2];
219 
220  deriv2[0] = 0.0;
221  if (setup.p >= 1)
222  deriv2[1] = 1.0;
223  for (int i=2; i<=setup.p; i++)
224  deriv2[i] = i*val2[i-1];
225  success = Stokhos::compareArrays(val1, "val1", val2, "val2",
226  setup.rtol, setup.atol, out);
227  success = success &&
228  Stokhos::compareArrays(deriv1, "deriv1", deriv2, "deriv2",
229  setup.rtol, setup.atol, out);
230 
231 
232  }
233 
234  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Evaluate ) {
235  double x = 0.1234;
236  Teuchos::Array<double> v1(setup.p+1), v2(setup.p+1);
237  for (int i=0; i<=setup.p; i++)
238  v1[i] = setup.basis.evaluate(x, i);
239 
240  // evaluate bases using formula
241  // He_0(x) = 1
242  // He_1(x) = x
243  // He_i(x) = x*He_{i-1}(x) - (i-1)*He_{i-2}(x), i=2,3,...
244  v2[0] = 1.0;
245  if (setup.p >= 1)
246  v2[1] = x;
247  for (int i=2; i<=setup.p; i++)
248  v2[i] = x*v2[i-1] - (i-1.0)*v2[i-2];
249  success = Stokhos::compareArrays(v1, "v1", v2, "v2", setup.rtol, setup.atol,
250  out);
251  }
252 
253  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, Recurrence ) {
254  Teuchos::Array<double> a1(setup.p+1), b1(setup.p+1), c1(setup.p+1),
255  d1(setup.p+1);
256  Teuchos::Array<double> a2(setup.p+1), b2(setup.p+1), c2(setup.p+1),
257  d2(setup.p+1);
258  setup.basis.getRecurrenceCoefficients(a1, b1, c1, d1);
259 
260  // compute coefficients using formula
261  a2[0] = 0.0; b2[0] = 1.0; c2[0] = 1.0; d2[0] = 1.0;
262  for (int i=1; i<=setup.p; i++) {
263  a2[i] = 0.0;
264  b2[i] = i;
265  c2[i] = 1.0;
266  d2[i] = 1.0;
267  }
268  success = true;
269  success = success &&
270  Stokhos::compareArrays(a1, "a1", a2, "a2", setup.rtol, setup.atol, out);
271  success = success &&
272  Stokhos::compareArrays(b1, "b1", b2, "b2", setup.rtol, setup.atol, out);
273  success = success &&
274  Stokhos::compareArrays(c1, "c1", c2, "c2", setup.rtol, setup.atol, out);
275  success = success &&
276  Stokhos::compareArrays(d1, "d1", d2, "d2", setup.rtol, setup.atol, out);
277  }
278 
279  // Tests alpha coefficients satisfy
280  // alpha_k = delta_k * (t*psi_k,psi_k)/(psi_k,psi_k)
281  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, RecurrenceAlpha ) {
282  Teuchos::Array<double> a1(setup.p+1), b1(setup.p+1), c1(setup.p+1),
283  d1(setup.p+1);
284  setup.basis.getRecurrenceCoefficients(a1, b1, c1, d1);
285 
287  const Teuchos::Array<double>& n1 = setup.basis.norm_squared();
290  setup.basis.getQuadPoints(2*setup.p, x, w, v);
291  int nqp = w.size();
292  for (int i=0; i<=setup.p; i++) {
293  a2[i] = 0;
294  for (int j=0; j<nqp; j++)
295  a2[i] += w[j]*x[j]*v[j][i]*v[j][i];
296  a2[i] = a2[i]*c1[i]/n1[i];
297  }
298  success = Stokhos::compareArrays(a1, "a1", a2, "a2", setup.rtol, setup.atol,
299  out);
300  }
301 
302  // Tests beta coefficients satisfy
303  // beta_k =
304  // gamma_k * delta_k/delta_{k-1} * (psi_k,psi_k)/(psi_{k-1},psi_{k-1})
305  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, RecurrenceBeta ) {
306  Teuchos::Array<double> a1(setup.p+1), b1(setup.p+1), c1(setup.p+1),
307  d1(setup.p+1);
308  setup.basis.getRecurrenceCoefficients(a1, b1, c1, d1);
309 
311  const Teuchos::Array<double>& n1 = setup.basis.norm_squared();
312  b2[0] = b1[0];
313  for (int i=1; i<=setup.p; i++) {
314  b2[i] = d1[i]*(c1[i]/c1[i-1])*(n1[i]/n1[i-1]);
315  }
316  success = Stokhos::compareArrays(b1, "b1", b2, "b2", setup.rtol, setup.atol,
317  out);
318  }
319 
320 #ifdef HAVE_STOKHOS_DAKOTA
321  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, QuadPointsDakota ) {
322  int n = static_cast<int>(std::ceil((setup.p+1)/2.0));
323  Teuchos::Array<double> x1, w1;
325  setup.basis.getQuadPoints(setup.p, x1, w1, v1);
326 
327  Teuchos::Array<double> x2(n), w2(n);
329  webbur::hermite_compute(n, &x2[0], &w2[0]);
330 
331  for (int i=0; i<n; i++) {
332  w2[i] *= 0.5/std::sqrt(std::atan(1.0)); // 1/sqrt(pi)
333  x2[i] *= std::sqrt(2.0);
334  v2[i].resize(setup.p+1);
335  setup.basis.evaluateBases(x2[i], v2[i]);
336  }
337  success = true;
338  success = success &&
339  Stokhos::compareArrays(x1, "x1", x2, "x2", setup.rtol, setup.atol, out);
340  success = success &&
341  Stokhos::compareArrays(w1, "w1", w2, "w2", setup.rtol, setup.atol, out);
342  for (int i=0; i<n; i++) {
343  std::stringstream ss1, ss2;
344  ss1 << "v1[" << i << "]";
345  ss2 << "v2[" << i << "]";
346  success = success &&
347  Stokhos::compareArrays(v1[i], ss1.str(), v2[i], ss2.str(),
348  setup.rtol, setup.atol, out);
349  }
350  }
351 #endif
352 
353 #ifdef HAVE_STOKHOS_FORUQTK
354  TEUCHOS_UNIT_TEST( Stokhos_HermiteBasis, QuadPointsForUQTK ) {
355  int n = static_cast<int>(std::ceil((setup.p+1)/2.0));
356  Teuchos::Array<double> x1, w1;
358  setup.basis.getQuadPoints(setup.p, x1, w1, v1);
359 
360  Teuchos::Array<double> x2(n), w2(n);
362  int kind = 4;
363  int kpts = 0;
364  double endpts[2] = {0.0, 0.0};
366  double alpha = 0.0;
367  double beta = 0.0;
368  GAUSSQ_F77(&kind, &n, &alpha, &beta, &kpts, endpts, &b[0], &x2[0], &w2[0]);
369 
370  for (int i=0; i<n; i++) {
371  w2[i] *= 0.5/std::sqrt(std::atan(1.0)); // 1/sqrt(pi)
372  x2[i] *= std::sqrt(2.0);
373  v2[i].resize(setup.p+1);
374  setup.basis.evaluateBases(x2[i], v2[i]);
375  }
376  success = true;
377  success = success &&
378  Stokhos::compareArrays(x1, "x1", x2, "x2", setup.rtol, setup.atol, out);
379  success = success &&
380  Stokhos::compareArrays(w1, "w1", w2, "w2", setup.rtol, setup.atol, out);
381  for (int i=0; i<n; i++) {
382  std::stringstream ss1, ss2;
383  ss1 << "v1[" << i << "]";
384  ss2 << "v2[" << i << "]";
385  success = success &&
386  Stokhos::compareArrays(v1[i], ss1.str(), v2[i], ss2.str(),
387  setup.rtol, setup.atol, out);
388  }
389  }
390 #endif
391 
392 }
393 
394 int main( int argc, char* argv[] ) {
395  Teuchos::GlobalMPISession mpiSession(&argc, &argv);
397 }
KOKKOS_INLINE_FUNCTION PCE< Storage > sqrt(const PCE< Storage > &a)
SerialBandDenseMatrixPrinter< OrdinalType, ScalarType > printMat(const SerialBandDenseMatrix< OrdinalType, ScalarType > &obj)
bool compareArrays(const Array1 &a1, const std::string &a1_name, const Array2 &a2, const std::string &a2_name, const ValueType &rel_tol, const ValueType &abs_tol, Teuchos::FancyOStream &out)
bool compareValues(const ValueType &a1, const std::string &a1_name, const ValueType &a2, const std::string &a2_name, const ValueType &rel_tol, const ValueType &abs_tol, Teuchos::FancyOStream &out)
static int runUnitTestsFromMain(int argc, char *argv[])
KOKKOS_INLINE_FUNCTION PCE< Storage > ceil(const PCE< Storage > &a)
UnitTestSetup< int, double > setup
void resize(size_type new_size, const value_type &x=value_type())
KOKKOS_INLINE_FUNCTION PCE< Storage > atan(const PCE< Storage > &a)
int main(int argc, char **argv)
ScalarTraits< ScalarType >::magnitudeType normInf() const
expr val()
#define TEUCHOS_TEST_EQUALITY(v1, v2, out, success)
Stokhos::HermiteBasis< OrdinalType, ValueType > basis
size_type size() const
TEUCHOS_UNIT_TEST(Stokhos_HermiteBasis, Order)
int n
void GAUSSQ_F77(int *kind, int *n, double *alpha, double *beta, int *kpts, double *endpts, double *b, double *x, double *w)