NeoPZ
Public Types | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
TPZGaussRule Class Reference

Implements the Gaussian quadrature. Numerical Integration Abstract class. More...

#include <tpzgaussrule.h>

Collaboration diagram for TPZGaussRule:
[legend]

Public Types

enum  { NRULESLOBATTO_ORDER, NRULESLEGENDRE_ORDER }
 

Public Member Functions

int NInt () const
 Returns number of integration points. More...
 
int Type ()
 Returns the gaussian quadrature type (Legendre, Lobatto, Jacobi) More...
 
long double Loc (int i) const
 Returns location of the ith point. More...
 
long double W (int i) const
 Returns weight for the ith point. More...
 
long double Alpha ()
 Returns ALPHA: it is the exponent of (1-X) in the quadrature rule: $ weight = (1-X)^\alpha * (1+X)^\beta $. More...
 
long double Beta ()
 Returns BETA: it is the exponent of (1-X) in the quadrature rule: $ weight = (1-X)^ALPHA * (1+X)^BETA $. More...
 
int Order ()
 return the order of the integration rule More...
 
void SetParametersJacobi (long double alpha, long double beta)
 Sets alpha and beta, the parameters of the Jacobi polynomials. More...
 
void SetType (int &type, int order)
 Sets a gaussian quadrature: 0 - Gauss Legendre, 1 - Gauss Lobatto, 2 - Gauss Jacobi with alpha = beta = 1., 3 - Gauss Chebyshev for alpha = beta = -0.5. More...
 
void Print (std::ostream &out=std::cout)
 Prints the number of integration points, all points and weights (as one dimension) More...
 

Protected Member Functions

int ComputingGaussLegendreQuadrature (int order)
 Computes the points and weights for Gauss Legendre Quadrature over the parametric 1D element [-1.0, 1.0]. More...
 
void ComputingGaussLegendreQuadrature (int *npoints, TPZVec< long double > &Location, TPZVec< long double > &Weight)
 Computes the points and weights for Gauss Legendre Quadrature over the parametric 1D element [-1.0, 1.0] It is util to be called by another rule that need the Gaussian quadrature based on number of points not on order. More...
 
int ComputingGaussLobattoQuadrature (int order)
 Computes the points and weights for Gauss Lobatto Quadrature over the parametric 1D element $ [-1.0, 1.0] $ It is to integrate functions $ f(x) $, but the first and last integration points are $ -1.0 $ and $1.0$ respectively. More...
 
int ComputingGaussJacobiQuadrature (int order, long double alpha, long double beta)
 Computes the points and weights for Gauss Lobatto Quadrature over the parametric 1D element $ [-1.0, 1.0] $ It is to integrate functions $ f(x) $, but one of integration point is $-1.0$ or $1.0$. More...
 
void ComputingGaussJacobiQuadrature (int *npoints, long double alpha, long double beta, TPZVec< long double > &Location, TPZVec< long double > &Weight)
 Compute numerical integration points and weight for Gauss Jacobi quadrature over $ [-1,1] $ It is to integrate functions as $ (1-x)^\alpha (1+x)^\beta f(x) $. It is util to be called by another rule. Do not stores the integration points and weights in the current instance. More...
 
long double JacobiPolinomial (long double x, int order, long double alpha, long double beta, TPZVec< long double > &b, TPZVec< long double > &c, long double *dp2, long double *p1)
 Evaluates the Jacobi polinomial for real x. More...
 
long double JacobiPolinomial (long double x, int alpha, int beta, unsigned int n)
 Evaluates the polinomial Jacobi on x point. More...
 
int ComputingGaussChebyshevQuadrature (int order)
 Computes the points and weights for Gauss Chebyshev Quadrature over the parametric 1D element [-1.0, 1.0] It is to integrate functions as $ f(x)/\sqrt(1-x^2)) $. More...
 

Private Member Functions

 TPZGaussRule (int order, int type=0, long double alpha=0.0L, long double beta=0.0L)
 Constructor of quadrature rule. More...
 
 TPZGaussRule (const TPZGaussRule &copy)
 
TPZGaussRuleoperator= (const TPZGaussRule &copy)
 
 ~TPZGaussRule ()
 Default destructor. More...
 

Private Attributes

int fType
 fType = 1 (Gauss Lobatto quadrature), fType = 2 (Gauss Jacobi quadrature) for any alpha and beta the parameters of Jacobi polynomials, fType = 3 (Gauss Chebyshev) is a special Gauss Jacobi quadrature for alpha = -.5 and beta = -0.5 fType = 0(default) (Gauss Legendre quadrature) is a special Gauss Jacobi quadrature for alpha = beta = 0.0 More...
 
int fNumInt
 Number of integration points for this object. More...
 
TPZVec< long double > fLocation
 Location of the integration point. More...
 
TPZVec< long double > fWeight
 Weight of the integration point. More...
 
long double fAlpha
 ALPHA is the exponent of (1-X) in the quadrature rule: weight = (1-X)^ALPHA * (1+X)^BETA. More...
 
long double fBeta
 BETA is the exponent of (1+X) in the quadrature rule: weight = (1-X)^ALPHA * (1+X)^BETA. More...
 
int fOrder
 order of the integration rule More...
 

Friends

class TPZIntRuleP3D
 A instance of the cubature rule for pyramid can access computing methods to construct its points and weights. More...
 
class TPZIntRuleList
 The list can to access the constructor of the current class. More...
 

Detailed Description

Implements the Gaussian quadrature. Numerical Integration Abstract class.

Author
Philippe R. B. Devloo phil@.nosp@m.fec..nosp@m.unica.nosp@m.mp.b.nosp@m.r

Definition at line 19 of file tpzgaussrule.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
NRULESLOBATTO_ORDER 
NRULESLEGENDRE_ORDER 

Definition at line 67 of file tpzgaussrule.h.

Constructor & Destructor Documentation

◆ TPZGaussRule() [1/2]

TPZGaussRule::TPZGaussRule ( int  order,
int  type = 0,
long double  alpha = 0.0L,
long double  beta = 0.0L 
)
private

Constructor of quadrature rule.

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule
typeIdentifies the type of integration points constructed
alphaExponent of the factor (1 - ksi) to Jacobi type
betaExponent of the factor (1 + ksi) to Jacobi type
Note
If alpha = beta = 0.0L the Jacobi type is identical to Legendre type

Definition at line 22 of file tpzgaussrule.cpp.

References IsZero(), and PZError.

◆ TPZGaussRule() [2/2]

TPZGaussRule::TPZGaussRule ( const TPZGaussRule copy)
private

Definition at line 66 of file tpzgaussrule.cpp.

◆ ~TPZGaussRule()

TPZGaussRule::~TPZGaussRule ( )
private

Default destructor.

Definition at line 88 of file tpzgaussrule.cpp.

References fLocation, fNumInt, fWeight, and TPZVec< T >::Resize().

Member Function Documentation

◆ Alpha()

long double TPZGaussRule::Alpha ( )
inline

Returns ALPHA: it is the exponent of (1-X) in the quadrature rule: $ weight = (1-X)^\alpha * (1+X)^\beta $.

Definition at line 81 of file tpzgaussrule.h.

References fAlpha.

◆ Beta()

long double TPZGaussRule::Beta ( )
inline

Returns BETA: it is the exponent of (1-X) in the quadrature rule: $ weight = (1-X)^ALPHA * (1+X)^BETA $.

Definition at line 83 of file tpzgaussrule.h.

References fBeta.

◆ ComputingGaussChebyshevQuadrature()

int TPZGaussRule::ComputingGaussChebyshevQuadrature ( int  order)
protected

Computes the points and weights for Gauss Chebyshev Quadrature over the parametric 1D element [-1.0, 1.0] It is to integrate functions as $ f(x)/\sqrt(1-x^2)) $.

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule

Definition at line 242 of file tpzgaussrule.cpp.

References fLocation, fNumInt, fWeight, and TPZVec< T >::Resize().

Referenced by SetParametersJacobi(), and SetType().

◆ ComputingGaussJacobiQuadrature() [1/2]

int TPZGaussRule::ComputingGaussJacobiQuadrature ( int  order,
long double  alpha,
long double  beta 
)
protected

Computes the points and weights for Gauss Lobatto Quadrature over the parametric 1D element $ [-1.0, 1.0] $ It is to integrate functions $ f(x) $, but one of integration point is $-1.0$ or $1.0$.

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule Compute numerical integration points and weight for Gauss Jacobi quadrature over $ [-1,1] $ It is to integrate functions as $ (1-x)^\alpha (1+x)^\beta f(x) $
orderOrder of the polinomial will be integrated exactly with this cubature rule
alphaExponent of the factor (1 - ksi) to Jacobi type
betaExponent of the factor (1 + ksi) to Jacobi type

Gauss Jacobi quadrature

Definition at line 379 of file tpzgaussrule.cpp.

References fLocation, fNumInt, and fWeight.

Referenced by TPZIntRuleP3D::ComputingCubatureRuleForPyramid(), SetParametersJacobi(), and SetType().

◆ ComputingGaussJacobiQuadrature() [2/2]

void TPZGaussRule::ComputingGaussJacobiQuadrature ( int *  npoints,
long double  alpha,
long double  beta,
TPZVec< long double > &  Location,
TPZVec< long double > &  Weight 
)
protected

Compute numerical integration points and weight for Gauss Jacobi quadrature over $ [-1,1] $ It is to integrate functions as $ (1-x)^\alpha (1+x)^\beta f(x) $. It is util to be called by another rule. Do not stores the integration points and weights in the current instance.

Parameters
npointsNumber of integration points
alphaExponent of factor (1+x)
betaExponent of factor (1-x)
LocationVector of integration points ksi
WeightVector of corresponding weights

Gauss Jacobi quadrature

Definition at line 390 of file tpzgaussrule.cpp.

References gamma(), JacobiPolinomial(), machinePrecision(), pow(), and TPZVec< T >::Resize().

◆ ComputingGaussLegendreQuadrature() [1/2]

int TPZGaussRule::ComputingGaussLegendreQuadrature ( int  order)
protected

Computes the points and weights for Gauss Legendre Quadrature over the parametric 1D element [-1.0, 1.0].

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule

GAUSS LEGENDRE QUADRATURE

Definition at line 170 of file tpzgaussrule.cpp.

References fLocation, fNumInt, fWeight, and NRULESLEGENDRE_ORDER.

Referenced by TPZIntRuleP3D::ComputingCubatureRuleForPyramid(), SetParametersJacobi(), and SetType().

◆ ComputingGaussLegendreQuadrature() [2/2]

void TPZGaussRule::ComputingGaussLegendreQuadrature ( int *  npoints,
TPZVec< long double > &  Location,
TPZVec< long double > &  Weight 
)
protected

Computes the points and weights for Gauss Legendre Quadrature over the parametric 1D element [-1.0, 1.0] It is util to be called by another rule that need the Gaussian quadrature based on number of points not on order.

Parameters
npointsNumber of integration points
LocationVector of integration points ksi
WeightVector of corresponding weights

Definition at line 183 of file tpzgaussrule.cpp.

References fabs, m, machinePrecision(), PZINTEGRAL_MAXITERATIONS_ALLOWED, TPZVec< T >::Resize(), and pzgeom::tol.

◆ ComputingGaussLobattoQuadrature()

int TPZGaussRule::ComputingGaussLobattoQuadrature ( int  order)
protected

Computes the points and weights for Gauss Lobatto Quadrature over the parametric 1D element $ [-1.0, 1.0] $ It is to integrate functions $ f(x) $, but the first and last integration points are $ -1.0 $ and $1.0$ respectively.

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule
Note
The implementation follow Karniadakis and Sherwin: "Spectral/hp element methods for computational fluid dynamics (Oxford University Press, 2005).

Computes the points and weights for Gauss Lobatto quadrature

Definition at line 288 of file tpzgaussrule.cpp.

References test::f, fabs, fLocation, fNumInt, fWeight, gamma(), JacobiPolinomial(), m, NRULESLOBATTO_ORDER, PZINTEGRAL_MAXITERATIONS_ALLOWED, TPZVec< T >::Resize(), and pzgeom::tol.

Referenced by SetParametersJacobi(), and SetType().

◆ JacobiPolinomial() [1/2]

long double TPZGaussRule::JacobiPolinomial ( long double  x,
int  order,
long double  alpha,
long double  beta,
TPZVec< long double > &  b,
TPZVec< long double > &  c,
long double *  dp2,
long double *  p1 
)
protected

Evaluates the Jacobi polinomial for real x.

Parameters
xReal number to compute J(x)
orderOrder of the polinomial
alphaExponent of factor (1+x)
betaExponent of factor (1-x)
bRecursion coefficients (Jacobi matrix)
cRecursion coefficients
dp2Returns the value of J'(x)
p1Returns the value of J[order-1](x), where J[order-1] represent the Jacobi polinomial to degree order-1.

Evaluate the Jacobi polinomial with alpha and beta real numbers

Definition at line 527 of file tpzgaussrule.cpp.

Referenced by ComputingGaussJacobiQuadrature(), ComputingGaussLobattoQuadrature(), and SetParametersJacobi().

◆ JacobiPolinomial() [2/2]

long double TPZGaussRule::JacobiPolinomial ( long double  x,
int  alpha,
int  beta,
unsigned int  n 
)
protected

Evaluates the polinomial Jacobi on x point.

Parameters
xpoint on the polinomial Jacobi is evaluated, J(x)
nis the order of the polinomial
alphais the exponent of the factor (1-x)
betais the exponent of the factor (1+x)
Note
The Jacobi polinomial is evaluated using a recursion formulae

Definition at line 262 of file tpzgaussrule.cpp.

◆ Loc()

long double TPZGaussRule::Loc ( int  i) const

Returns location of the ith point.

Definition at line 96 of file tpzgaussrule.cpp.

References fLocation, fNumInt, and PZError.

Referenced by TPZPrInteg< TFather >::Point(), and Type().

◆ NInt()

int TPZGaussRule::NInt ( ) const
inline

Returns number of integration points.

Definition at line 70 of file tpzgaussrule.h.

References fNumInt.

Referenced by TPZPrInteg< TFather >::NPoints(), and TPZPrInteg< TFather >::Point().

◆ operator=()

TPZGaussRule & TPZGaussRule::operator= ( const TPZGaussRule copy)
private

Definition at line 72 of file tpzgaussrule.cpp.

References fAlpha, fBeta, fLocation, fNumInt, fOrder, fType, and fWeight.

◆ Order()

int TPZGaussRule::Order ( )
inline

return the order of the integration rule

Definition at line 86 of file tpzgaussrule.h.

References fOrder.

Referenced by TPZInt1d::SetOrder(), TPZIntCube3D::SetOrder(), TPZInt1d::TPZInt1d(), TPZIntCube3D::TPZIntCube3D(), and TPZIntQuad::TPZIntQuad().

◆ Print()

void TPZGaussRule::Print ( std::ostream &  out = std::cout)

Prints the number of integration points, all points and weights (as one dimension)

Definition at line 116 of file tpzgaussrule.cpp.

References fLocation, fNumInt, fType, and fWeight.

Referenced by TPZInt1d::Print(), and SetParametersJacobi().

◆ SetParametersJacobi()

void TPZGaussRule::SetParametersJacobi ( long double  alpha,
long double  beta 
)
inline

Sets alpha and beta, the parameters of the Jacobi polynomials.

Parameters
alphaExponent of the factor (1 - ksi) to Jacobi type
betaExponent of the factor (1 + ksi) to Jacobi type

Definition at line 95 of file tpzgaussrule.h.

References ComputingGaussChebyshevQuadrature(), ComputingGaussJacobiQuadrature(), ComputingGaussLegendreQuadrature(), ComputingGaussLobattoQuadrature(), gamma(), JacobiPolinomial(), machinePrecision(), Print(), and SetType().

◆ SetType()

void TPZGaussRule::SetType ( int &  type,
int  order 
)

Sets a gaussian quadrature: 0 - Gauss Legendre, 1 - Gauss Lobatto, 2 - Gauss Jacobi with alpha = beta = 1., 3 - Gauss Chebyshev for alpha = beta = -0.5.

Parameters
orderOrder of the polinomial will be integrated exactly with this cubature rule
typeIdentifies the type of integration points constructed

Definition at line 131 of file tpzgaussrule.cpp.

References ComputingGaussChebyshevQuadrature(), ComputingGaussJacobiQuadrature(), ComputingGaussLegendreQuadrature(), ComputingGaussLobattoQuadrature(), DebugStop, fAlpha, fBeta, fLocation, fType, fWeight, and TPZVec< T >::Resize().

Referenced by SetParametersJacobi(), TPZInt1d::SetType(), TPZIntQuad::SetType(), and TPZIntCube3D::SetType().

◆ Type()

int TPZGaussRule::Type ( )
inline

Returns the gaussian quadrature type (Legendre, Lobatto, Jacobi)

Definition at line 73 of file tpzgaussrule.h.

References fType, Loc(), and W().

Referenced by TPZIntQuad::SetOrder().

◆ W()

long double TPZGaussRule::W ( int  i) const

Returns weight for the ith point.

Definition at line 107 of file tpzgaussrule.cpp.

References fNumInt, fWeight, and PZError.

Referenced by TPZPrInteg< TFather >::Point(), and Type().

Friends And Related Function Documentation

◆ TPZIntRuleList

friend class TPZIntRuleList
friend

The list can to access the constructor of the current class.

Definition at line 32 of file tpzgaussrule.h.

◆ TPZIntRuleP3D

friend class TPZIntRuleP3D
friend

A instance of the cubature rule for pyramid can access computing methods to construct its points and weights.

Definition at line 29 of file tpzgaussrule.h.

Member Data Documentation

◆ fAlpha

long double TPZGaussRule::fAlpha
private

ALPHA is the exponent of (1-X) in the quadrature rule: weight = (1-X)^ALPHA * (1+X)^BETA.

Definition at line 42 of file tpzgaussrule.h.

Referenced by Alpha(), operator=(), and SetType().

◆ fBeta

long double TPZGaussRule::fBeta
private

BETA is the exponent of (1+X) in the quadrature rule: weight = (1-X)^ALPHA * (1+X)^BETA.

Definition at line 44 of file tpzgaussrule.h.

Referenced by Beta(), operator=(), and SetType().

◆ fLocation

TPZVec<long double> TPZGaussRule::fLocation
private

◆ fNumInt

int TPZGaussRule::fNumInt
private

◆ fOrder

int TPZGaussRule::fOrder
private

order of the integration rule

Definition at line 47 of file tpzgaussrule.h.

Referenced by operator=(), and Order().

◆ fType

int TPZGaussRule::fType
private

fType = 1 (Gauss Lobatto quadrature), fType = 2 (Gauss Jacobi quadrature) for any alpha and beta the parameters of Jacobi polynomials, fType = 3 (Gauss Chebyshev) is a special Gauss Jacobi quadrature for alpha = -.5 and beta = -0.5 fType = 0(default) (Gauss Legendre quadrature) is a special Gauss Jacobi quadrature for alpha = beta = 0.0

Definition at line 26 of file tpzgaussrule.h.

Referenced by operator=(), Print(), SetType(), and Type().

◆ fWeight

TPZVec<long double> TPZGaussRule::fWeight
private

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