|
IPPL (Independent Parallel Particle Layer)
IPPL
|
This is class represents the Gauss-Legendre quadrature rule. It is a special case of the Gauss-Jacobi quadrature rule with alpha = beta = 0.0. More...
#include <GaussJacobiQuadrature.h>
Public Types | |
| using | scalar_t = double |
Public Member Functions | |
| 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. | |
| void | computeNodesAndWeights () override |
| scalar_t | getChebyshevNodes (const size_t &i) const |
| Returns the i-th Chebyshev node, used as initial guess for the Newton iterations. | |
| size_t | getOrder () const |
| Returns the order of the quadrature rule. (order = degree + 1). | |
| size_t | getDegree () const |
| Returns the degree of exactness of the quadrature rule. | |
| Vector< T, numElementNodes > | getWeightsForRefElement () const |
| Get the quadrature weights for the reference element. | |
| Vector< Vector< T, dim >, numElementNodes > | getIntegrationNodesForRefElement () const |
| Get the integration (quadrature) nodes for the reference element. | |
| Vector< T, NumNodes1D > | getIntegrationNodes1D (const T &a, const T &b) const |
| Get the quadrature nodes for one dimension. (With respect to the given domain [a, b]). | |
| Vector< T, NumNodes1D > | getWeights1D (const T &a, const T &b) const |
| Get the quadrature weights for one dimension. (With respect to the given domain [a, b]). | |
Static Public Attributes | |
| static constexpr unsigned | numNodes1D = NumNodes1D |
| static constexpr unsigned | dim = ElementType::dim |
| static constexpr unsigned | numElementNodes |
Protected Attributes | |
| unsigned | degree_m |
| const ElementType & | ref_element_m |
| Vector< T, NumNodes1D > | integration_nodes_m |
| Vector< T, NumNodes1D > | weights_m |
| T | a_m |
| T | b_m |
Private Member Functions | |
| 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 from LehrFEM++. | |
Private Attributes | |
| const T | alpha_m |
| const T | beta_m |
| const size_t | max_newton_iterations_m |
| const size_t | min_newton_iterations_m |
This is class represents the Gauss-Legendre quadrature rule. It is a special case of the Gauss-Jacobi quadrature rule with alpha = beta = 0.0.
| T | floating point number type of the quadrature nodes and weights |
| NumNodes1D | number of quadrature nodes for one dimension |
| ElementType | element type for which the quadrature rule is defined |
Definition at line 90 of file GaussJacobiQuadrature.h.
|
inherited |
Definition at line 31 of file GaussJacobiQuadrature.h.
|
inline |
Construct a new Gauss Legendre Quadrature rule object.
| ref_element | reference element to compute the quadrature nodes on |
| max_newton_itersations | maximum number of Newton iterations (default 10) |
| min_newton_iterations | minimum number of Newton iterations (default 1) |
Definition at line 99 of file GaussJacobiQuadrature.h.
References ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature().
|
overridevirtualinherited |
Computes the quadrature nodes and weights and stores them in the quadrature nodes and weights arrays.
the following algorithm for computing the roots and weights for the Gauss-Jacobi quadrature has been taken from LehrFEM++ (MIT License) https://craffael.github.io/lehrfempp/gauss__quadrature_8cc_source.html
Implements ippl::Quadrature< T, NumNodes1D, ElementType >.
Definition at line 91 of file GaussJacobiQuadrature.hpp.
References alpha_m, beta_m, ippl::Chebyshev, getChebyshevNodes(), getLehrFEMInitialGuess(), ippl::Quadrature< T, NumNodes1D, ElementType >::integration_nodes_m, ippl::LehrFEM, max_newton_iterations_m, min_newton_iterations_m, and ippl::Quadrature< T, NumNodes1D, ElementType >::weights_m.
Referenced by GaussJacobiQuadrature().
|
inherited |
Returns the i-th Chebyshev node, used as initial guess for the Newton iterations.
| i | index of the Chebyshev node |
Definition at line 34 of file GaussJacobiQuadrature.hpp.
Referenced by computeNodesAndWeights().
|
inherited |
Returns the degree of exactness of the quadrature rule.
Definition at line 13 of file Quadrature.hpp.
References degree_m.
|
inherited |
Get the quadrature nodes for one dimension. (With respect to the given domain [a, b]).
| a | local domain start |
| b | local domain end |
Definition at line 72 of file Quadrature.hpp.
Referenced by getIntegrationNodesForRefElement().
|
inherited |
Get the integration (quadrature) nodes for the reference element.
Definition at line 47 of file Quadrature.hpp.
References getIntegrationNodes1D(), and numElementNodes.
|
privateinherited |
Computes the initial guess for the Newton iterations, the way they are computed in the implementation from LehrFEM++.
| i | index of the initial guess (corresponding to the i-th quadrature node) |
| integration_nodes | the integration nodes |
Definition at line 41 of file GaussJacobiQuadrature.hpp.
References alpha_m, and beta_m.
Referenced by computeNodesAndWeights().
|
inherited |
Returns the order of the quadrature rule. (order = degree + 1).
Definition at line 8 of file Quadrature.hpp.
References degree_m.
|
inherited |
Get the quadrature weights for one dimension. (With respect to the given domain [a, b]).
| a | local domain start |
| b | local domain end |
Definition at line 81 of file Quadrature.hpp.
References a_m, b_m, and weights_m.
Referenced by getWeightsForRefElement().
|
inherited |
Get the quadrature weights for the reference element.
Definition at line 19 of file Quadrature.hpp.
References getWeights1D(), and numElementNodes.
|
protectedinherited |
Definition at line 117 of file Quadrature.h.
Referenced by ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature(), getWeights1D(), and ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::MidpointQuadrature().
|
privateinherited |
Definition at line 74 of file GaussJacobiQuadrature.h.
Referenced by computeNodesAndWeights(), GaussJacobiQuadrature(), and getLehrFEMInitialGuess().
|
protectedinherited |
Definition at line 119 of file Quadrature.h.
Referenced by ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature(), getWeights1D(), and ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::MidpointQuadrature().
|
privateinherited |
Definition at line 75 of file GaussJacobiQuadrature.h.
Referenced by computeNodesAndWeights(), GaussJacobiQuadrature(), and getLehrFEMInitialGuess().
|
protectedinherited |
Definition at line 111 of file Quadrature.h.
Referenced by ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature(), getDegree(), getOrder(), and ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::MidpointQuadrature().
|
staticconstexprinherited |
Definition at line 41 of file Quadrature.h.
|
protectedinherited |
Definition at line 113 of file Quadrature.h.
Referenced by ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), and ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature().
|
privateinherited |
Definition at line 77 of file GaussJacobiQuadrature.h.
Referenced by computeNodesAndWeights(), and GaussJacobiQuadrature().
|
privateinherited |
Definition at line 78 of file GaussJacobiQuadrature.h.
Referenced by computeNodesAndWeights(), and GaussJacobiQuadrature().
|
staticconstexprinherited |
Definition at line 44 of file Quadrature.h.
Referenced by getIntegrationNodesForRefElement(), and getWeightsForRefElement().
|
staticconstexprinherited |
Definition at line 38 of file Quadrature.h.
|
protectedinherited |
Definition at line 112 of file Quadrature.h.
Referenced by Quadrature().
|
protectedinherited |
Definition at line 114 of file Quadrature.h.
Referenced by ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), ippl::MidpointQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights(), ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature(), and getWeights1D().