-
Notifications
You must be signed in to change notification settings - Fork 10
/
Monegato.h
380 lines (335 loc) · 13.6 KB
/
Monegato.h
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
#ifndef EIGEN_MONEGATO_H
#define EIGEN_MONEGATO_H
namespace Eigen
{
/**
* \ingroup NumericalIntegration_Module
*
* \class Monegato
*
* \brief A class for computing the Gauss/Kronrod weights.
*
* \tparam Scalar floating point type
*
* This class is based on "Some remarks on the construction of
* extended Gaussian quadrature rules", Giovanni Monegato,
* Math. Comp., Vol. 32 (1978) pp. 247-252. http://www.jstor.org/stable/2006272 .
*
* The code is based on quadpackcpp (http://quadpackpp.sourceforge.net/)
* library written by Jerry Gagelman <[email protected]> which is distributed under
* GNU General Public License
*/
template <typename Scalar>
class Monegato
{
public:
/**
* \brief Compute the absolute value of a scalar
* \param[in] x inpute scalar value
* \return absolute value of scalar input
*/
static Scalar abs(Scalar x)
{
return (x < Scalar(0)) ? -(x) : x;
}
/**
* \brief Compute the Legendre polynomials and their error
* \param[in] n degree
* \param[in] x value of the variable
* \param[out] err error to be returned
* \return value of the Legendre polynomial
*
* This function is based on the routine gsl_sf_legendre_Pl_e
* distributed with GSL
*/
static Scalar legendre_err(const int n, const Scalar x, Scalar& err)
{
if (n == 0)
{
err = Scalar(0);
return Scalar(1);
}
else if (n == 1)
{
err = Scalar(0);
return x;
}
// below Sree modified this to avoid -Wmaybe-uninitialized
Scalar P0 = Scalar(1), P1 = x, P2=x;
Scalar E0 = NumTraits<Scalar>::epsilon();
Scalar E1 = abs(x) * NumTraits<Scalar>::epsilon();
for (int k = 1; k < n; ++k)
{
P2 = ((2*k + 1) * x * P1 - k * P0) / (k + 1);
err = ((2*k + 1) * abs(x) * E1 + k * E0) / (2*(k + 1));
P0 = P1; P1 = P2;
E0 = E1; E1 = err;
}
return P2;
}
/**
* \brief Three-term recursion identity for the Legendre derivatives
* \param[in] n degree
* \param[in] x value of the variable
* \return value of the Legendre derivative
*/
static Scalar legendre_deriv(int const n, Scalar const x)
{
if (n == 0)
{
return Scalar(0);
}
else if (n == 1)
{
return Scalar(1);
}
Scalar P0 = Scalar(1);
Scalar P1 = x;
Scalar P2 = P1;
Scalar dP0 = Scalar(0);
Scalar dP1 = Scalar(1);
Scalar dP2 = P1;
for (int k = 1; k < n; ++k)
{
P2 = ((2*k + 1) * x * P1 - k * P0) / (k + Scalar(1));
dP2 = (2*k + 1) * P1 + dP0;
P0 = P1;
P1 = P2;
dP0 = dP1;
dP1 = dP2;
}
return dP2;
}
/**
* \brief Compute derivatives of Chebyshev polynomials
* \param[in] x variable
* \param[in] n_ degree
* \param[in] coefs Chebyshev coefficients
* \return Derivative of Chebyshev polynomials
*/
static Scalar chebyshev_series_deriv(const Scalar x,
const int n_,
const Eigen::Array<Scalar,
Eigen::Dynamic, 1>& coefs)
{
Scalar d1(0), d2(0);
Scalar y2 = 2 * x; // linear term for Clenshaw recursion
for (int k = n_; k >= 2; --k)
{
Scalar temp = d1;
d1 = y2 * d1 - d2 + k * coefs(k);
d2 = temp;
}
return y2 * d1 - d2 + coefs(1);
}
/**
* \breif Evaluation of the Chebyshev polynomial using Clenshaw recursion
* \param[in] x variable
* \param[in] n_ degree
* \param[in] coefs Chebyshev coefficients
* \param[out] err error to be returned
* \return value of Chebyshev polynomial
*/
static Scalar chebyshev_series(const Scalar x,
const int n_,
const Eigen::Array<Scalar, Eigen::Dynamic, 1>& coefs,
Scalar& err)
{
Scalar d1(0), d2(0);
Scalar absc = abs(coefs(0)); // final term for truncation error
Scalar y2 = 2 * x; // linear term for Clenshaw recursion
for (int k = n_; k >= 1; --k)
{
Scalar temp = d1;
d1 = y2 * d1 - d2 + coefs(k);
d2 = temp;
absc += abs(coefs(k));
}
err = absc * NumTraits<Scalar>::epsilon();
return x * d1 - d2 + coefs(0)/2.;
}
/**
* \brief Computes the zeros of the Legendre polynomial
* \param[in] m_ degree
* \param[out] zeros the zeros of the Legendre polynomial
*/
static void legendre_zeros(const int m_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& zeros)
{
Eigen::Array<Scalar, Eigen::Dynamic, 1> temp = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(m_+1);
zeros(0) = Scalar(-1);
zeros(1) = Scalar(1);
Scalar delta, epsilon;
for (int k = 1; k <= m_; ++k)
{
// Loop to locate zeros of P_k interlacing z_0,...,z_k
for (int j = 0; j < k; ++j)
{
// Newton's method for P_k :
// initialize solver at midpoint of (z_j, z_{j+1})
delta = 1;
Scalar x_j = (zeros(j) + zeros(j+1)) / 2.;
Scalar P_k = legendre_err(k, x_j, epsilon);
while (abs(P_k) > epsilon &&
abs(delta) > NumTraits<Scalar>::epsilon())
{
delta = P_k / legendre_deriv(k, x_j);
x_j -= delta;
P_k = legendre_err(k, x_j, epsilon);
}
temp(j) = x_j;
}
// Copy roots tmp_0,...,tmp_{k-1} to z_1,...z_k:
zeros(k+1) = zeros(k);
for (int j = 0; j < k; ++j)
{
zeros(j+1) = temp(j);
}
}
}
/**
* \brief Computes coefficients of the Chebyshev polynomial.
* \param[in] m_ degree of the Chebyshev polynomial.
* \param[out] coefs coefficients of the Chebyshev polynomial.
*/
static void chebyshev_coefs(const int m_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& coefs)
{
int ell = (m_ + 1)/2;
Eigen::Array<Scalar, Eigen::Dynamic, 1> alpha = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(ell+1);
Eigen::Array<Scalar, Eigen::Dynamic, 1> f = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(ell+1);
// Care must be exercised in initalizing the constants in the definitions.
// Compilers interpret expressions like "(2*k + 1.0)/(k + 2.0)" as floating
// point precision, before casting to Real.
f(1) = Scalar(m_+1) / Scalar(2*m_ + 3);
alpha(0) = Scalar(1); // coefficient of T_{m+1}
alpha(1) = -f(1);
for (int k = 2; k <= ell; ++k)
{
f(k) = f(k-1) * (2*k - 1) * (m_ + k) / (k * (2*m_ + 2*k + 1));
alpha(k) = -f(k);
for (int i = 1; i < k; ++i)
{
alpha(k) -= f(i) * alpha(k-i);
}
}
for (int k = 0; k <= ell; ++k)
{
coefs(m_ + 1 - 2*k) = alpha(k);
if (m_ >= 2*k)
{
coefs(m_ - 2*k) = Scalar(0);
}
}
}
/**
* \brief Compute the Gauss/Kronrod abscissae
* \param[in] n_ size of Gauss-Kronrod arrays
* \param[in] Gauss-Legendre degree
* \param[in] coefs Chebyshev coefficients
* \param[out] xgk_ Gauss/Kronrod abscissae
*/
static void gauss_kronrod_abscissae(const int n_,
const int m_,
const Eigen::Array<Scalar, Eigen::Dynamic, 1>& zeros,
const Eigen::Array<Scalar, Eigen::Dynamic, 1>& coefs,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& xgk_)
{
// Now from the function gauss_kronrod_abscissae
Scalar delta, epsilon;
for (int k = 0; k < n_ / 2; ++k)
{
delta = 1;
// Newton's method for E_{n+1} :
Scalar x_k = (zeros(m_-k) + zeros(m_+1-k))/Scalar(2);
Scalar E = chebyshev_series(x_k,n_,coefs, epsilon);
while (abs(E) > epsilon &&
abs(delta) > NumTraits<Scalar>::epsilon() )
{
delta = E / chebyshev_series_deriv(x_k,n_,coefs);
x_k -= delta;
E = chebyshev_series(x_k,n_,coefs, epsilon);
}
xgk_(2*k) = x_k;
// Copy adjacent Legendre-zero into the array:
if (2*k+1 < n_)
{
xgk_(2*k+1) = zeros(m_-k);
}
}
}
/**
* \brief Compute Gauss/Kronrod weights
* \param[in] n_ size of Gauss-Kronrod arrays
* \param[in] m_ Gauss-Legendre degree
* \param[in] xgk_ Gauss/Kronrod abscissae
* \param[out] wg_ Gauss weights
* \param[out] wgk_ Kronrod weights
*/
static void gauss_kronrod_weights(const int& n_,
const int m_,
const Eigen::Array<Scalar, Eigen::Dynamic, 1>& coefs,
const Eigen::Array<Scalar, Eigen::Dynamic, 1>& xgk_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& wg_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& wgk_)
{
Scalar err;
// Gauss-Legendre weights:
for(int k = 0; k < n_ / 2; ++k)
{
Scalar x = xgk_(2*k + 1);
wg_(k) = (Scalar(-2) / ((m_ + 1) * legendre_deriv(m_, x) * legendre_err(m_+1, x, err)));
}
// The ratio of leading coefficients of P_n and T_{n+1} is computed
// from the recursive formulae for the respective polynomials.
Scalar F_m = Scalar(2) / Scalar(2*m_ + 1);
for (int k = 1; k <= m_; ++k)
{
F_m *= (Scalar(2*k) / Scalar(2*k - 1));
}
// Gauss-Kronrod weights:
for (int k = 0; k < n_; ++k)
{
Scalar x = xgk_(k);
if (k % 2 == 0)
{
wgk_(k) = F_m / (legendre_err(m_, x, err) * chebyshev_series_deriv(x,n_,coefs));
}
else
{
wgk_(k) = (wg_(k/2) + F_m / (legendre_deriv(m_, x) * chebyshev_series(x,n_,coefs, err)));
}
}
}
/**
* \brief A method ofr computing Gauss/Kronrod abscissae and weights
* \param[in] m_ Gauss-Legendre degree
* \param[out] xgk_ Gauss/Kronrod abscissae
* \param[out] wgk_ Gauss/Kronrod weights
* \param[out] xk_ Gauss abscissae
* \param[out] wg_ Gauss weights
*
* Note that Gauss abscissae is not calculated here.
*/
static void computeAbscissaeAndWeights(unsigned int m_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& xgk_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& wgk_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& xk_,
Eigen::Array<Scalar, Eigen::Dynamic, 1>& wg_)
{
const unsigned int n_ = m_ + 1;
xgk_ = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(n_); //2*nNodes+1
wgk_ = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(n_); //2*nNodes+1
xk_ = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(n_/2); //2*nNodes
wg_ = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(n_/2); //2*nNodes
// initialise the coefficients to zero
Eigen::Array<Scalar, Eigen::Dynamic, 1> coefs = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(n_+1);
Eigen::Array<Scalar, Eigen::Dynamic, 1> zeros = Eigen::Array<Scalar, Eigen::Dynamic, 1>::Zero(m_+2);
legendre_zeros(m_, zeros);
chebyshev_coefs(m_, coefs);
gauss_kronrod_abscissae(n_, m_, zeros,coefs, xgk_);
gauss_kronrod_weights(n_, m_,coefs, xgk_, wg_, wgk_ );
}
};
} // namespace Eigen
#endif // EIGEN_MONEGATO_H