10#include "gauss_quadrature.h"
12#include <boost/math/special_functions/gamma.hpp>
13#include <boost/multiprecision/cpp_bin_float.hpp>
17 unsigned num_points) {
18 LF_ASSERT_MSG(num_points > 0,
"num_points must be positive.");
21 boost::multiprecision::number<boost::multiprecision::cpp_bin_float<
22 57, boost::multiprecision::digit_base_2>>;
24 static const scalar_t kPi = boost::math::constants::pi<scalar_t>();
26 Eigen::VectorXd points(num_points);
27 Eigen::VectorXd weights(num_points);
31 unsigned int m = (num_points + 1) / 2;
34 for (
unsigned int i = 0; i < m; ++i) {
36 scalar_t z = cos(kPi * (i + 0.75) / (num_points + 0.5));
47 for (
unsigned int j = 0; j < num_points; ++j) {
50 p1 = ((2.0 * j + 1.) * z * p2 - j * p3) / (j + 1.);
56 pp = num_points * (z * p1 - p2) / (z * z - 1.0);
60 }
while (abs(z - z1) > 1e-17);
62 points(i) = (0.5 * (1 - z)).convert_to<
double>();
63 points(num_points - 1 - i) = (0.5 * (1 + z)).convert_to<
double>();
64 weights(i) = (1. / ((1.0 - z * z) * pp * pp)).convert_to<double>();
65 weights(num_points - 1 - i) = weights(i);
68 return {std::move(points), std::move(weights)};
73 LF_ASSERT_MSG(num_points > 0,
"num_points must be positive.");
74 LF_ASSERT_MSG(alpha > -1,
"alpha > -1 required");
75 LF_ASSERT_MSG(beta > -1,
"beta > -1 required.");
77 using boost::math::lgamma;
79 boost::multiprecision::number<boost::multiprecision::cpp_bin_float<
80 57, boost::multiprecision::digit_base_2>>;
101 std::vector<scalar_t> points(num_points);
102 Eigen::VectorXd weights(num_points);
108 an = alpha / num_points;
109 bn = beta / num_points;
111 (2.78 / (4.0 + num_points * num_points) + 0.768 * an / num_points);
112 r2 = 1.0 + 1.48 * an + 0.96 * bn + 0.452 * an * an + 0.83 * an * bn;
116 r1 = (4.1 + alpha) / ((1.0 + alpha) * (1.0 + 0.156 * alpha));
117 r2 = 1.0 + 0.06 * (num_points - 8.0) * (1.0 + 0.12 * alpha) / num_points;
118 r3 = 1.0 + 0.012 * beta * (1.0 + 0.25 * std::abs(alpha)) / num_points;
119 z -= (1.0 - z) * r1 * r2 * r3;
122 r1 = (1.67 + 0.28 * alpha) / (1.0 + 0.37 * alpha);
123 r2 = 1.0 + 0.22 * (num_points - 8.0) / num_points;
124 r3 = 1.0 + 8.0 * beta / ((6.28 + beta) * num_points * num_points);
125 z -= (points[0] - z) * r1 * r2 * r3;
126 }
else if (i == num_points - 2) {
128 r1 = (1.0 + 0.235 * beta) / (0.766 + 0.119 * beta);
129 r2 = 1.0 / (1.0 + 0.639 * (num_points - 4.0) /
130 (1.0 + 0.71 * (num_points - 4.0)));
132 (1.0 + 20.0 * alpha / ((7.5 + alpha) * num_points * num_points));
133 z += (z - points[num_points - 4]) * r1 * r2 * r3;
134 }
else if (i == num_points - 1) {
136 r1 = (1.0 + 0.37 * beta) / (1.67 + 0.28 * beta);
137 r2 = 1.0 / (1.0 + 0.22 * (num_points - 8.0) / num_points);
139 (1.0 + 8.0 * alpha / ((6.28 + alpha) * num_points * num_points));
140 z += (z - points[num_points - 3]) * r1 * r2 * r3;
143 z = 3.0 * points[i - 1] - 3.0 * points[i - 2] + points[i - 3];
145 alfbet = alpha + beta;
147 for (its = 1; its <= MAX_IT; ++its) {
153 p1 = (alpha - beta + temp * z) / 2.0;
158 temp = 2 * j + alfbet;
159 a = 2 * j * (j + alfbet) * (temp - 2.0);
161 (alpha * alpha - beta * beta + temp * (temp - 2.0) * z);
162 c = 2.0 * (j - 1 + alpha) * (j - 1 + beta) * temp;
163 p1 = (b * p2 - c * p3) / a;
165 pp = (num_points * (alpha - beta - temp * z) * p1 +
166 2.0 * (num_points + alpha) * (num_points + beta) * p2) /
167 (temp * (1.0 - z * z));
173 if (abs(z - z1) <= 1e-17) {
177 LF_VERIFY_MSG(its <= MAX_IT,
"too many iterations.");
180 weights(i) = (exp(lgamma(alpha + num_points) + lgamma(beta + num_points) -
181 lgamma(num_points + 1.) -
182 lgamma(
double(num_points) + alfbet + 1.0)) *
183 temp * pow(2.0, alfbet) / (pp * p2))
184 .convert_to<
double>();
187 Eigen::VectorXd points_result(num_points);
189 points_result(i) = points[i].convert_to<
double>();
192 return {points_result, weights};
Rules for numerical quadrature on reference entity shapes.
std::tuple< Eigen::VectorXd, Eigen::VectorXd > GaussJacobi(quadDegree_t num_points, double alpha, double beta)
Computes the quadrature points and weights for the interval [-1,1] of a Gauss-Jacobi quadrature rule ...
std::tuple< Eigen::VectorXd, Eigen::VectorXd > GaussLegendre(unsigned num_points)
unsigned int quadDegree_t