IPPL (Independent Parallel Particle Layer)
IPPL
Loading...
Searching...
No Matches
ippl::GaussLegendreQuadrature< T, NumNodes1D, ElementType > Class Template Reference

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>

Inheritance diagram for ippl::GaussLegendreQuadrature< T, NumNodes1D, ElementType >:
Collaboration diagram for ippl::GaussLegendreQuadrature< T, NumNodes1D, ElementType >:

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, numElementNodesgetWeightsForRefElement () const
 Get the quadrature weights for the reference element.
Vector< Vector< T, dim >, numElementNodesgetIntegrationNodesForRefElement () 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

Detailed Description

template<typename T, unsigned NumNodes1D, typename ElementType>
class ippl::GaussLegendreQuadrature< T, NumNodes1D, ElementType >

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.

Template Parameters
Tfloating point number type of the quadrature nodes and weights
NumNodes1Dnumber of quadrature nodes for one dimension
ElementTypeelement type for which the quadrature rule is defined

Definition at line 90 of file GaussJacobiQuadrature.h.

Member Typedef Documentation

◆ scalar_t

template<typename T, unsigned NumNodes1D, typename ElementType>
using ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::scalar_t = double
inherited

Definition at line 31 of file GaussJacobiQuadrature.h.

Constructor & Destructor Documentation

◆ GaussLegendreQuadrature()

template<typename T, unsigned NumNodes1D, typename ElementType>
ippl::GaussLegendreQuadrature< T, NumNodes1D, ElementType >::GaussLegendreQuadrature ( const ElementType & ref_element,
const size_t & max_newton_itersations = 10,
const size_t & min_newton_iterations = 1 )
inline

Construct a new Gauss Legendre Quadrature rule object.

Parameters
ref_elementreference element to compute the quadrature nodes on
max_newton_itersationsmaximum number of Newton iterations (default 10)
min_newton_iterationsminimum number of Newton iterations (default 1)

Definition at line 99 of file GaussJacobiQuadrature.h.

References ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::GaussJacobiQuadrature().

Here is the call graph for this function:

Member Function Documentation

◆ computeNodesAndWeights()

template<typename T, unsigned NumNodes1D, typename ElementType>
void ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::computeNodesAndWeights ( )
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().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getChebyshevNodes()

template<typename T, unsigned NumNodes1D, typename ElementType>
GaussJacobiQuadrature< T, NumNodes1D, ElementType >::scalar_t ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::getChebyshevNodes ( const size_t & i) const
inherited

Returns the i-th Chebyshev node, used as initial guess for the Newton iterations.

Parameters
iindex of the Chebyshev node
Returns
scalar_t - i-th Chebyshev node

Definition at line 34 of file GaussJacobiQuadrature.hpp.

Referenced by computeNodesAndWeights().

Here is the caller graph for this function:

◆ getDegree()

template<typename T, unsigned NumNodes1D, typename ElementType>
size_t ippl::Quadrature< T, NumNodes1D, ElementType >::getDegree ( ) const
inherited

Returns the degree of exactness of the quadrature rule.

Returns
unsigned - degree

Definition at line 13 of file Quadrature.hpp.

References degree_m.

◆ getIntegrationNodes1D()

template<typename T, unsigned NumNodes1D, typename ElementType>
Vector< T, NumNodes1D > ippl::Quadrature< T, NumNodes1D, ElementType >::getIntegrationNodes1D ( const T & a,
const T & b ) const
inherited

Get the quadrature nodes for one dimension. (With respect to the given domain [a, b]).

Parameters
alocal domain start
blocal domain end
Returns
Vector<T, NumNodes1D>

Definition at line 72 of file Quadrature.hpp.

Referenced by getIntegrationNodesForRefElement().

Here is the caller graph for this function:

◆ getIntegrationNodesForRefElement()

template<typename T, unsigned NumNodes1D, typename ElementType>
Vector< Vector< T, Quadrature< T, NumNodes1D, ElementType >::dim >, Quadrature< T, NumNodes1D, ElementType >::numElementNodes > ippl::Quadrature< T, NumNodes1D, ElementType >::getIntegrationNodesForRefElement ( ) const
inherited

Get the integration (quadrature) nodes for the reference element.

Returns
Vector<Vector<T, Dim>, numElementNodes>

Definition at line 47 of file Quadrature.hpp.

References getIntegrationNodes1D(), and numElementNodes.

Here is the call graph for this function:

◆ getLehrFEMInitialGuess()

template<typename T, unsigned NumNodes1D, typename ElementType>
GaussJacobiQuadrature< T, NumNodes1D, ElementType >::scalar_t ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::getLehrFEMInitialGuess ( const size_t & i,
const Vector< scalar_t, NumNodes1D > & integration_nodes ) const
privateinherited

Computes the initial guess for the Newton iterations, the way they are computed in the implementation from LehrFEM++.

Parameters
iindex of the initial guess (corresponding to the i-th quadrature node)
integration_nodesthe integration nodes
Returns
scalar_t - initial guess

Definition at line 41 of file GaussJacobiQuadrature.hpp.

References alpha_m, and beta_m.

Referenced by computeNodesAndWeights().

Here is the caller graph for this function:

◆ getOrder()

template<typename T, unsigned NumNodes1D, typename ElementType>
size_t ippl::Quadrature< T, NumNodes1D, ElementType >::getOrder ( ) const
inherited

Returns the order of the quadrature rule. (order = degree + 1).

Returns
unsigned - order

Definition at line 8 of file Quadrature.hpp.

References degree_m.

◆ getWeights1D()

template<typename T, unsigned NumNodes1D, typename ElementType>
Vector< T, NumNodes1D > ippl::Quadrature< T, NumNodes1D, ElementType >::getWeights1D ( const T & a,
const T & b ) const
inherited

Get the quadrature weights for one dimension. (With respect to the given domain [a, b]).

Parameters
alocal domain start
blocal domain end
Returns
Vector<T, NumNodes1D>

Definition at line 81 of file Quadrature.hpp.

References a_m, b_m, and weights_m.

Referenced by getWeightsForRefElement().

Here is the caller graph for this function:

◆ getWeightsForRefElement()

template<typename T, unsigned NumNodes1D, typename ElementType>
Vector< T, Quadrature< T, NumNodes1D, ElementType >::numElementNodes > ippl::Quadrature< T, NumNodes1D, ElementType >::getWeightsForRefElement ( ) const
inherited

Get the quadrature weights for the reference element.

Returns
Vector<T, numElementNodes>

Definition at line 19 of file Quadrature.hpp.

References getWeights1D(), and numElementNodes.

Here is the call graph for this function:

Member Data Documentation

◆ a_m

◆ alpha_m

template<typename T, unsigned NumNodes1D, typename ElementType>
const T ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::alpha_m
privateinherited

◆ b_m

◆ beta_m

template<typename T, unsigned NumNodes1D, typename ElementType>
const T ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::beta_m
privateinherited

◆ degree_m

template<typename T, unsigned NumNodes1D, typename ElementType>
unsigned ippl::Quadrature< T, NumNodes1D, ElementType >::degree_m
protectedinherited

◆ dim

template<typename T, unsigned NumNodes1D, typename ElementType>
unsigned ippl::Quadrature< T, NumNodes1D, ElementType >::dim = ElementType::dim
staticconstexprinherited

Definition at line 41 of file Quadrature.h.

◆ integration_nodes_m

template<typename T, unsigned NumNodes1D, typename ElementType>
Vector<T, NumNodes1D> ippl::Quadrature< T, NumNodes1D, ElementType >::integration_nodes_m
protectedinherited

◆ max_newton_iterations_m

template<typename T, unsigned NumNodes1D, typename ElementType>
const size_t ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::max_newton_iterations_m
privateinherited

Definition at line 77 of file GaussJacobiQuadrature.h.

Referenced by computeNodesAndWeights(), and GaussJacobiQuadrature().

◆ min_newton_iterations_m

template<typename T, unsigned NumNodes1D, typename ElementType>
const size_t ippl::GaussJacobiQuadrature< T, NumNodes1D, ElementType >::min_newton_iterations_m
privateinherited

Definition at line 78 of file GaussJacobiQuadrature.h.

Referenced by computeNodesAndWeights(), and GaussJacobiQuadrature().

◆ numElementNodes

template<typename T, unsigned NumNodes1D, typename ElementType>
unsigned ippl::Quadrature< T, NumNodes1D, ElementType >::numElementNodes
staticconstexprinherited
Initial value:

Definition at line 44 of file Quadrature.h.

Referenced by getIntegrationNodesForRefElement(), and getWeightsForRefElement().

◆ numNodes1D

template<typename T, unsigned NumNodes1D, typename ElementType>
unsigned ippl::Quadrature< T, NumNodes1D, ElementType >::numNodes1D = NumNodes1D
staticconstexprinherited

Definition at line 38 of file Quadrature.h.

◆ ref_element_m

template<typename T, unsigned NumNodes1D, typename ElementType>
const ElementType& ippl::Quadrature< T, NumNodes1D, ElementType >::ref_element_m
protectedinherited

Definition at line 112 of file Quadrature.h.

Referenced by Quadrature().

◆ weights_m


The documentation for this class was generated from the following file: