IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
GaussJacobiQuadrature.h
Go to the documentation of this file.
1// Class GaussJacobiQuadrature
2// The GaussJacobiQuadrature class. This is a class representing a Gauss-Jacobi quadrature
3// The algoirthm in computeNodesAndWeights is based on the LehrFEM++ library.
4// https://craffael.github.io/lehrfempp/gauss__quadrature_8cc_source.html
5
6#ifndef IPPL_GAUSSJACOBIQUADRATURE_H
7#define IPPL_GAUSSJACOBIQUADRATURE_H
8
10
11namespace ippl {
12
17
26 template <typename T, unsigned NumNodes1D, typename ElementType>
27 class GaussJacobiQuadrature : public Quadrature<T, NumNodes1D, ElementType> {
28 public:
29 // a higher precision floating point number type used for the computation
30 // of the quadrature nodes and weights
31 using scalar_t = double; // might be equivalant to double, depending on compiler
32
42 GaussJacobiQuadrature(const ElementType& ref_element, const T& alpha, const T& beta,
43 const size_t& max_newton_itersations = 10,
44 const size_t& min_newton_iterations = 1);
45
50 void computeNodesAndWeights() override;
51
59 scalar_t getChebyshevNodes(const size_t& i) const; // FIXME maybe move somewhere else?
60
61 private:
72 const size_t& i, const Vector<scalar_t, NumNodes1D>& integration_nodes) const;
73
74 const T alpha_m;
75 const T beta_m;
76
79 };
80
89 template <typename T, unsigned NumNodes1D, typename ElementType>
90 class GaussLegendreQuadrature : public GaussJacobiQuadrature<T, NumNodes1D, ElementType> {
91 public:
99 GaussLegendreQuadrature(const ElementType& ref_element,
100 const size_t& max_newton_itersations = 10,
101 const size_t& min_newton_iterations = 1)
102 : GaussJacobiQuadrature<T, NumNodes1D, ElementType>(
103 ref_element, 0.0, 0.0, max_newton_itersations, min_newton_iterations) {}
104 };
105
114 template <typename T, unsigned NumNodes1D, typename ElementType>
115 class ChebyshevGaussQuadrature : public GaussJacobiQuadrature<T, NumNodes1D, ElementType> {
116 public:
124 ChebyshevGaussQuadrature(const ElementType& ref_element,
125 const size_t& max_newton_itersations = 10,
126 const size_t& min_newton_iterations = 1)
127 : GaussJacobiQuadrature<T, NumNodes1D, ElementType>(
128 ref_element, -0.5, -0.5, max_newton_itersations, min_newton_iterations) {}
129 };
130
131} // namespace ippl
132
134
135#endif
Definition Archive.h:20
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.
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...
GaussLegendreQuadrature(const ElementType &ref_element, const size_t &max_newton_itersations=10, const size_t &min_newton_iterations=1)
Construct a new Gauss Legendre Quadrature rule object.
ChebyshevGaussQuadrature(const ElementType &ref_element, const size_t &max_newton_itersations=10, const size_t &min_newton_iterations=1)
Construct a new Chebyshev Gauss Quadrature rule object.
Quadrature(const ElementType &ref_element)
Construct a new Quadrature object.
Definition Quadrature.hpp:4