5 template <
typename T,
unsigned NumNodes1D,
typename ElementType>
8 const size_t& max_newton_iterations,
const size_t& min_newton_iterations)
12 , max_newton_iterations_m(max_newton_iterations)
13 , min_newton_iterations_m(min_newton_iterations) {
14 assert(
alpha > -1.0 &&
"alpha >= -1.0 is not satisfied");
15 assert(beta > -1.0 &&
"beta >= -1.0 is not satisfied");
16 assert(max_newton_iterations >= 1 &&
"max_newton_iterations >= 1 is not satisfied");
19 &&
"min_newton_iterations_m <= max_newton_iterations_m is not satisfied");
32 template <
typename T,
unsigned NumNodes1D,
typename ElementType>
35 return -Kokkos::cos((2.0 *
static_cast<scalar_t>(i) + 1.0) * Kokkos::numbers::pi_v<scalar_t>
36 / (2.0 * NumNodes1D));
39 template <
typename T,
unsigned NumNodes1D,
typename ElementType>
44 integration_nodes)
const {
50 scalar_t z = (i > 0) ? integration_nodes[i - 1] : 0.0;
55 const scalar_t bn = beta / NumNodes1D;
56 r1 = (1.0 +
alpha) * (2.78 / (4.0 + NumNodes1D * NumNodes1D) + 0.768 * an / NumNodes1D);
57 r2 = 1.0 + 1.48 * an + 0.96 * bn + 0.452 * an * an + 0.83 * an * bn;
62 r2 = 1.0 + 0.06 * (NumNodes1D - 8.0) * (1.0 + 0.12 *
alpha) / NumNodes1D;
63 r3 = 1.0 + 0.012 * beta * (1.0 + 0.25 * Kokkos::abs(
alpha)) / NumNodes1D;
64 z -= (1.0 - z) * r1 * r2 * r3;
67 r1 = (1.67 + 0.28 *
alpha) / (1.0 + 0.37 *
alpha);
68 r2 = 1.0 + 0.22 * (NumNodes1D - 8.0) / NumNodes1D;
69 r3 = 1.0 + 8.0 * beta / ((6.28 + beta) * NumNodes1D * NumNodes1D);
70 z -= (integration_nodes[0] - z) * r1 * r2 * r3;
71 }
else if (i == NumNodes1D - 2) {
73 r1 = (1.0 + 0.235 * beta) / (0.766 + 0.119 * beta);
74 r2 = 1.0 / (1.0 + 0.639 * (NumNodes1D - 4.0) / (1.0 + 0.71 * (NumNodes1D - 4.0)));
75 r3 = 1.0 / (1.0 + 20.0 *
alpha / ((7.5 +
alpha) * NumNodes1D * NumNodes1D));
76 z += (z - integration_nodes[NumNodes1D - 4]) * r1 * r2 * r3;
77 }
else if (i == NumNodes1D - 1) {
79 r1 = (1.0 + 0.37 * beta) / (1.67 + 0.28 * beta);
80 r2 = 1.0 / (1.0 + 0.22 * (NumNodes1D - 8.0) / NumNodes1D);
81 r3 = 1.0 / (1.0 + 8.0 *
alpha / ((6.28 +
alpha) * NumNodes1D * NumNodes1D));
82 z += (z - integration_nodes[NumNodes1D - 3]) * r1 * r2 * r3;
85 z = 3.0 * z - 3.0 * integration_nodes[i - 2] + integration_nodes[i - 3];
90 template <
typename T,
unsigned NumNodes1D,
typename ElementType>
124 for (
size_t i = 0; i < NumNodes1D; ++i) {
127 z = this->getLehrFEMInitialGuess(i, integration_nodes);
129 z = -this->getChebyshevNodes(i);
131 throw IpplException(
"GaussJacobiQuadrature::computeNodesAndWeights",
132 "Unknown initial guess type");
145 p1 = (
alpha - beta + temp * z) / 2.0;
147 for (
size_t j = 2; j <= NumNodes1D; ++j) {
150 temp = 2 * j + alfbet;
151 a = 2 * j * (j + alfbet) * (temp - 2.0);
152 b = (temp - 1.0) * (
alpha *
alpha - beta * beta + temp * (temp - 2.0) * z);
153 c = 2.0 * (j - 1 +
alpha) * (j - 1 + beta) * temp;
154 p1 = (b * p2 -
c * p3) / a;
156 pp = (NumNodes1D * (
alpha - beta - temp * z) * p1
157 + 2.0 * (NumNodes1D +
alpha) * (NumNodes1D + beta) * p2)
158 / (temp * (1.0 - z * z));
168 if (its > this->min_newton_iterations_m && Kokkos::abs(z - z1) <= tolerance) {
172 }
while (its <= this->max_newton_iterations_m);
174 if (its > this->max_newton_iterations_m) {
176 std::cout <<
"Root " << NumNodes1D - i - 1
177 <<
" didn't converge. Tolerance may be too high for data type"
184 integration_nodes[i] = z;
188 Kokkos::exp(Kokkos::lgamma(
alpha + NumNodes1D) + Kokkos::lgamma(beta + NumNodes1D)
189 - Kokkos::lgamma(NumNodes1D + 1.)
190 - Kokkos::lgamma(
static_cast<double>(NumNodes1D) + alfbet + 1.0))
191 * temp * Kokkos::pow(2.0, alfbet) / (pp * p2);
194 this->integration_nodes_m[i] =
static_cast<T>(-integration_nodes[i]);
195 this->weights_m[i] =
static_cast<T>(weights[i]);
Inform & endl(Inform &inf)
constexpr double alpha
The fine structure constant, no dimension.
constexpr double e
The value of.
constexpr double c
The velocity of light in m/s.
Implementations for FFT constructor/destructor and transforms.
scalar_t getChebyshevNodes(const size_t &i) const
Returns the i-th Chebyshev node, used as initial guess for the Newton iterations.
GaussJacobiQuadrature(const ElementType &ref_element, const T &alpha, const T &beta, const size_t &max_newton_itersations=10, const size_t &min_newton_iterations=1)
Construct a new Gauss Jacobi Quadrature rule object.
const size_t min_newton_iterations_m
scalar_t getLehrFEMInitialGuess(const size_t &i, const Vector< scalar_t, NumNodes1D > &integration_nodes) const
Computes the initial guess for the Newton iterations, the way they are computed in the implementation...
const size_t max_newton_iterations_m
void computeNodesAndWeights() override
Computes the quadrature nodes and weights and stores them in the quadrature nodes and weights arrays.
This is the base class for all quadrature rules.
Vector< T, NumNodes1D > integration_nodes_m
Vector< T, NumNodes1D > weights_m