Nektar++
TestPolylib.cpp
Go to the documentation of this file.
1///////////////////////////////////////////////////////////////////////////////
2//
3// File: TestPolylib.cpp
4//
5// For more information, please see: http://www.nektar.info
6//
7// The MIT License
8//
9// Copyright (c) 2006 Division of Applied Mathematics, Brown University (USA),
10// Department of Aeronautics, Imperial College London (UK), and Scientific
11// Computing and Imaging Institute, University of Utah (USA).
12//
13// Permission is hereby granted, free of charge, to any person obtaining a
14// copy of this software and associated documentation files (the "Software"),
15// to deal in the Software without restriction, including without limitation
16// the rights to use, copy, modify, merge, publish, distribute, sublicense,
17// and/or sell copies of the Software, and to permit persons to whom the
18// Software is furnished to do so, subject to the following conditions:
19//
20// The above copyright notice and this permission notice shall be included
21// in all copies or substantial portions of the Software.
22//
23// THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
24// OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
25// FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL
26// THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
27// LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING
28// FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER
29// DEALINGS IN THE SOFTWARE.
30//
31// Description: Unit tests for the Polylib library.
32//
33///////////////////////////////////////////////////////////////////////////////
34
36
37#include <boost/test/tools/floating_point_comparison.hpp>
38#include <boost/test/unit_test.hpp>
39
40#include <iostream>
41#include <tuple>
42#include <vector>
43
45{
46
47/// Lower range of number of quadrature points to use for unit tests.
48#define NPLOWER 5
49/// Upper range of number of quadrature points to use for unit tests.
50#define NPUPPER 15
51/// Tolerance to be used for floating point comparisons.
52#define EPS 1e-11
53
54/// Dot product two vectors.
55double ddot(int n, double *x, int incx, double *y, int incy)
56{
57 double sum = 0.0;
58
59 while (n--)
60 {
61 sum += (*x) * (*y);
62 x += incx;
63 y += incy;
64 }
65
66 return sum;
67}
68
69// Storage
70std::vector<double> z(NPUPPER), w(NPUPPER), p(NPUPPER);
71std::vector<double> d(NPUPPER *NPUPPER), q(NPUPPER *NPUPPER);
72
73/**
74 * @brief Test integrals of Jacobi polynomials.
75 *
76 * This routine evaluates the integral \f[ \int_{-1}^1 (1-x)^\alpha (1+x)^\beta
77 * P_n^{\alpha,\beta} dx = 0 \f] using \f$ -0.5 \leq \alpha,\beta \leq 5 \f$ and
78 * using \f$ N \f$ points where \f$ N \f$ lies between #NPLOWER and #NPUPPER.
79 * Tolerance is taken using the #EPS value.
80 */
81inline void TestIntegral(
82 std::function<void(double *, double *, int, double, double)> func, int o)
83{
84 double alpha = -0.5;
85 while (alpha <= 5.0)
86 {
87 double beta = -0.5;
88 while (beta <= 5.0)
89 {
90 for (int np = NPLOWER; np <= NPUPPER; ++np)
91 {
92 func(&z[0], &w[0], np, alpha, beta);
93 for (int n = 2; n < 2 * np - 1 - o; ++n)
94 {
95 Polylib::jacobfd(np, &z[0], &p[0], nullptr, n, alpha, beta);
96 double sum = ddot(np, &w[0], 1, &p[0], 1);
97 BOOST_CHECK_SMALL(sum, EPS);
98 }
99 }
100
101 beta += 0.5;
102 }
103
104 alpha += 0.5;
105 }
106}
107
108/**
109 * @brief Test derivatives using Gaussian quadrature.
110 *
111 * This routine evaluates the deriatives \f[ \frac{d}{dx}(x^n) = nx^{n-1} \f]
112 * using \f$ -0.5 \leq \alpha,\beta \leq 5 \f$ and using \f$ N \f$ points where
113 * \f$ N \f$ lies between #NPLOWER and #NPUPPER. Tolerance is taken using the
114 * #EPS value.
115 */
117 std::function<void(double *, double *, int, double, double)> func,
118 std::function<void(double *, double *, int, double, double)> funcD)
119{
120 double alpha = -0.5;
121 while (alpha <= 5.0)
122 {
123 double beta = -0.5;
124 while (beta <= 5.0)
125 {
126 for (int np = NPLOWER; np <= NPUPPER; ++np)
127 {
128 func(&z[0], &w[0], np, alpha, beta);
129 funcD(&d[0], &z[0], np, alpha, beta);
130 for (int n = 2; n < np - 1; ++n)
131 {
132 for (int i = 0; i < np; ++i)
133 {
134 p[i] = pow(z[i], n);
135 }
136
137 double sum = 0.0;
138
139 for (int i = 0; i < np; ++i)
140 {
141 sum += fabs(ddot(np, &d[0] + i, np, &p[0], 1) -
142 n * pow(z[i], n - 1));
143 }
144 sum /= np;
145
146 BOOST_CHECK_SMALL(sum, EPS);
147 }
148 }
149 beta += 0.5;
150 }
151 alpha += 0.5;
152 }
153}
154
155/**
156 * @brief Evaluate interpolation using Gaussian quadrature.
157 *
158 * This routine evaluates the interpolation of \f$ z^n \f$ to \f$ x^n \f$, where
159 * \f$ z \f$ are the quadrature zeros and \f$ x \f$ are the equispaced points
160 * \f[ x_n = \frac{2n}{N-1} - 1, \quad 0\leq n\leq N, \f] using \f$ -0.5 \leq
161 * \alpha,\beta \leq 5 \f$ and using \f$ N \f$ points where \f$ N \f$ lies
162 * between #NPLOWER and #NPUPPER. Tolerance is taken using the
163 * #EPS value.
164 */
166 std::function<void(double *, double *, int, double, double)> func,
167 std::function<void(double *, double *, double *, int, int, double, double)>
168 funcI)
169{
170 double alpha = -0.5;
171 while (alpha <= 5.0)
172 {
173 double beta = -0.5;
174 while (beta <= 5.0)
175 {
176 for (int np = NPLOWER; np <= NPUPPER; ++np)
177 {
178 func(&z[0], &w[0], np, alpha, beta);
179 for (int n = 2; n < np - 1; ++n)
180 {
181 for (int i = 0; i < np; ++i)
182 {
183 w[i] = 2.0 * i / (double)(np - 1) - 1.0;
184 p[i] = pow(z[i], n);
185 }
186 funcI(&d[0], &z[0], &w[0], np, np, alpha, beta);
187
188 double sum = 0.0;
189 for (int i = 0; i < np; ++i)
190 {
191 sum += fabs(ddot(np, &d[0] + i, np, &p[0], 1) -
192 pow(w[i], n));
193 }
194 sum /= np;
195
196 BOOST_CHECK_SMALL(sum, EPS);
197 }
198 }
199 beta += 0.5;
200 }
201 alpha += 0.5;
202 }
203}
204
206 std::function<void(double *, double *, int, double, double)> func,
207 std::function<void(double *, double *, int)> funcQ)
208{
209 double alpha = -0.5;
210 while (alpha <= 5.0)
211 {
212 double beta = -0.5;
213 while (beta <= 5.0)
214 {
215 for (int np = NPLOWER; np <= NPUPPER; ++np)
216 {
217 func(&z[0], &w[0], np, alpha, beta);
218 funcQ(&q[0], &z[0], np);
219 for (int n = 2; n < np - 1; ++n)
220 {
221 for (int i = 0; i < np; ++i)
222 {
223 p[i] = pow(z[i], n);
224 }
225
226 double sum = 0.0;
227
228 for (int i = 0; i < np; ++i)
229 {
230 sum += fabs(ddot(np, &q[0] + i * np, 1, &p[0], 1) -
231 pow(z[i], n + 1) / (n + 1));
232 }
233 sum /= np;
234
235 BOOST_CHECK_SMALL(sum, EPS);
236 }
237 }
238 beta += 0.5;
239 }
240 alpha += 0.5;
241 }
242}
243
245{
247}
248
249BOOST_AUTO_TEST_CASE(TestGaussRadauM)
250{
252}
253
254BOOST_AUTO_TEST_CASE(TestGaussRadauP)
255{
257}
258
259BOOST_AUTO_TEST_CASE(TestGaussLobatto)
260{
262}
263
265{
267}
268
269BOOST_AUTO_TEST_CASE(TestGaussRadauMDiff)
270{
272}
273
274BOOST_AUTO_TEST_CASE(TestGaussRadauPDiff)
275{
277}
278
279BOOST_AUTO_TEST_CASE(TestGaussLobattoDiff)
280{
282}
283
284BOOST_AUTO_TEST_CASE(TestGaussInterp)
285{
287}
288
289BOOST_AUTO_TEST_CASE(TestGaussRadauMInterp)
290{
292}
293
294BOOST_AUTO_TEST_CASE(TestGaussRadauPInterp)
295{
297}
298
299BOOST_AUTO_TEST_CASE(TestGaussLobattoInterp)
300{
302}
303
304BOOST_AUTO_TEST_CASE(TestGaussIntMatrix)
305{
307}
308
309BOOST_AUTO_TEST_CASE(TestGaussRadauMIntMatrix)
310{
312}
313
314BOOST_AUTO_TEST_CASE(TestGaussRadauPIntMatrix)
315{
317}
318
319BOOST_AUTO_TEST_CASE(TestGaussLobattoIntMatrix)
320{
322}
323
324BOOST_AUTO_TEST_CASE(TestGammaFraction)
325{
326 double a, c;
327
328 using arrT = std::tuple<int, double, int, double>;
329 std::vector<std::tuple<double, arrT, arrT>> ac = {
330 {362880.0, std::make_tuple(10, 0.0, 0, 1.0),
331 std::make_tuple(1, 0.0, 12, -2.0)},
332 {30.0, std::make_tuple(4, 3., 5, 0.), std::make_tuple(5, 0., 7, 0.)},
333 {21651.09375, std::make_tuple(7, 3.5, 5, 0.5),
334 std::make_tuple(5, 0.5, 10, 0.5)},
335 {97429.921875, std::make_tuple(10, 0.5, 5, -0.5),
336 std::make_tuple(5, -0.5, 11, -0.5)},
337 {2279.0625, std::make_tuple(10, -0.5, 5, 0.5),
338 std::make_tuple(5, 0.5, 10, -0.5)},
339 {639383.8623046875, std::make_tuple(10, 0.5, 0, 0.5),
340 std::make_tuple(0, 0.5, 10, 0.5)},
341 {5967498288235982848., std::make_tuple(100, 0., 90, 0.5),
342 std::make_tuple(90, 0.5, 100, 0.)},
343 {5618641603777298169856., std::make_tuple(200, 0., 191, -0.5),
344 std::make_tuple(190, 0.5, 200, 0.)},
345 {77396694214720029196288., std::make_tuple(200, 0., 190, 0.),
346 std::make_tuple(190, 0., 200, 0.)}};
347
348 for (auto &test : ac)
349 {
350 double a = std::get<0>(test);
351 arrT c1 = std::get<1>(test);
352 arrT c2 = std::get<2>(test);
353 double c1e = Polylib::gammaFracGammaF(std::get<0>(c1), std::get<1>(c1),
354 std::get<2>(c1), std::get<3>(c1));
355 double c2e = Polylib::gammaFracGammaF(std::get<0>(c2), std::get<1>(c2),
356 std::get<2>(c2), std::get<3>(c2));
357
358 BOOST_CHECK_SMALL(c1e / a - 1.0, EPS);
359 BOOST_CHECK_SMALL(c2e * a - 1.0, EPS);
360 }
361
362 int Q = 2;
363 for (double alpha = -0.5; alpha <= 150.; alpha += 0.5)
364 {
365 for (double beta = -0.5; beta <= 150.; beta += 0.5)
366 {
367 a = Polylib::gammaF(Q + alpha) / Polylib::gammaF(Q + beta);
368 c = Polylib::gammaFracGammaF(Q, alpha, Q, beta) / a - 1.;
369 BOOST_CHECK_SMALL(c, EPS);
370 }
371 }
372}
373
374} // namespace Nektar::UnitTests
#define NPUPPER
Upper range of number of quadrature points to use for unit tests.
Definition: TestPolylib.cpp:50
#define EPS
Tolerance to be used for floating point comparisons.
Definition: TestPolylib.cpp:52
#define NPLOWER
Lower range of number of quadrature points to use for unit tests.
Definition: TestPolylib.cpp:48
@ beta
Gauss Radau pinned at x=-1,.
Definition: PointsType.h:59
std::vector< double > w(NPUPPER)
std::vector< double > p(NPUPPER)
std::vector< double > z(NPUPPER)
void TestInterpolation(std::function< void(double *, double *, int, double, double)> func, std::function< void(double *, double *, double *, int, int, double, double)> funcI)
Evaluate interpolation using Gaussian quadrature.
std::vector< double > d(NPUPPER *NPUPPER)
double ddot(int n, double *x, int incx, double *y, int incy)
Dot product two vectors.
Definition: TestPolylib.cpp:55
void TestIntegral(std::function< void(double *, double *, int, double, double)> func, int o)
Test integrals of Jacobi polynomials.
Definition: TestPolylib.cpp:81
std::vector< double > q(NPUPPER *NPUPPER)
void TestIntegralMatrix(std::function< void(double *, double *, int, double, double)> func, std::function< void(double *, double *, int)> funcQ)
void TestDifferentiation(std::function< void(double *, double *, int, double, double)> func, std::function< void(double *, double *, int, double, double)> funcD)
Test derivatives using Gaussian quadrature.
BOOST_AUTO_TEST_CASE(TestGaussInt)
double gammaF(const double x)
Calculate the Gamma function , , for integer values and halves.
Definition: Polylib.cpp:1413
void Dgj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.
Definition: Polylib.cpp:653
void zwgrjm(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=-1.
Definition: Polylib.cpp:182
void zwgrjp(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Radau-Jacobi zeros and weights with end point at z=1.
Definition: Polylib.cpp:225
void Dglj(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Lobatto-Jacobi zeros.
Definition: Polylib.cpp:833
void Imgrjp(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at ...
Definition: Polylib.cpp:1116
void zwglj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.
Definition: Polylib.cpp:266
void zwgj(double *z, double *w, const int np, const double alpha, const double beta)
Gauss-Jacobi zeros and weights.
Definition: Polylib.cpp:152
void Imgj(double *im, const double *zgj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Jacobi points to an arbitrary distribution at points zm.
Definition: Polylib.cpp:1054
void Qg(double *Q, const double *z, const int np)
Compute the Integration Matrix.
Definition: Polylib.cpp:611
double gammaFracGammaF(const int x, const double alpha, const int y, const double beta)
Calculate fraction of two Gamma functions, , for integer values and halves.
Definition: Polylib.cpp:1468
void Imgrjm(double *im, const double *zgrj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at...
Definition: Polylib.cpp:1085
void Dgrjm(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a z...
Definition: Polylib.cpp:704
void Dgrjp(double *D, const double *z, const int np, const double alpha, const double beta)
Compute the Derivative Matrix associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.
Definition: Polylib.cpp:768
void Imglj(double *im, const double *zglj, const double *zm, const int nz, const int mz, const double alpha, const double beta)
Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.
Definition: Polylib.cpp:1147
void jacobfd(const int np, const double *z, double *poly_in, double *polyd, const int n, const double alpha, const double beta)
Routine to calculate Jacobi polynomials, , and their first derivative, .
Definition: Polylib.cpp:1248