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