OGS
NumLib::IntegrationGaussLegendreTri Class Reference

Detailed Description

Gauss-Legendre quadrature rule for triangles.

Gauss-Legendre quadrature rule for triangles is originally given as

\[ \int F(x,y) dx dy = \int F(x(r, s), y(r, s)) j(r,s) dr ds \approx \frac{1}{2} \sum_i ( F(x(r_i, s_i), y(r_i, s_i)) w_i ) \]

To make it consistent with other elements, we rewrite the above formula as

\[ \int F(x,y) dx dy \approx \sum_i ( F(x(r_i, s_i), y(r_i, s_i)) w'_i ) \]

by defining the new weight \( w'=\frac{1}{2} w \).

Definition at line 36 of file IntegrationGaussLegendreTri.h.

#include <IntegrationGaussLegendreTri.h>

Public Member Functions

 IntegrationGaussLegendreTri (unsigned order=2)
 
void setIntegrationOrder (unsigned order)
 Change the integration order.
 
unsigned getIntegrationOrder () const
 return current integration order.
 
unsigned getNumberOfPoints () const
 return the number of sampling points
 
MathLib::WeightedPoint getWeightedPoint (unsigned const igp) const
 

Static Public Member Functions

static MathLib::WeightedPoint getWeightedPoint (unsigned const order, unsigned const igp)
 
template<typename Method >
static MathLib::WeightedPoint getWeightedPoint (unsigned const igp)
 
static unsigned getNumberOfPoints (unsigned order)
 

Private Attributes

unsigned _order
 
unsigned _n_sampl_pt {0}
 

Constructor & Destructor Documentation

◆ IntegrationGaussLegendreTri()

NumLib::IntegrationGaussLegendreTri::IntegrationGaussLegendreTri ( unsigned order = 2)
inlineexplicit

Construct this object with the given integration order

Parameters
orderintegration order (default 2)

Definition at line 44 of file IntegrationGaussLegendreTri.h.

44 : _order(order)
45 {
46 this->setIntegrationOrder(order);
47 }
void setIntegrationOrder(unsigned order)
Change the integration order.

References setIntegrationOrder().

Member Function Documentation

◆ getIntegrationOrder()

unsigned NumLib::IntegrationGaussLegendreTri::getIntegrationOrder ( ) const
inline

return current integration order.

Definition at line 57 of file IntegrationGaussLegendreTri.h.

57{ return _order; }

References _order.

Referenced by getWeightedPoint().

◆ getNumberOfPoints() [1/2]

unsigned NumLib::IntegrationGaussLegendreTri::getNumberOfPoints ( ) const
inline

return the number of sampling points

Definition at line 60 of file IntegrationGaussLegendreTri.h.

References _n_sampl_pt.

Referenced by setIntegrationOrder().

◆ getNumberOfPoints() [2/2]

static unsigned NumLib::IntegrationGaussLegendreTri::getNumberOfPoints ( unsigned order)
inlinestatic

get the number of integration points

Parameters
orderthe number of integration points
Returns
the number of points

Definition at line 109 of file IntegrationGaussLegendreTri.h.

110 {
111 switch (order)
112 {
113 case 1:
115 case 2:
117 case 3:
119 case 4:
121 }
122 OGS_FATAL("Integration order {:d} not implemented for triangles.",
123 order);
124 }
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

◆ getWeightedPoint() [1/3]

template<typename Method >
static MathLib::WeightedPoint NumLib::IntegrationGaussLegendreTri::getWeightedPoint ( unsigned const igp)
inlinestatic

Definition at line 98 of file IntegrationGaussLegendreTri.h.

99 {
100 return MathLib::WeightedPoint(Method::X[igp], 0.5 * Method::W[igp]);
101 }

◆ getWeightedPoint() [2/3]

MathLib::WeightedPoint NumLib::IntegrationGaussLegendreTri::getWeightedPoint ( unsigned const igp) const
inline

get coordinates of a integration point

Parameters
igpThe integration point index
Returns
a weighted point

Definition at line 68 of file IntegrationGaussLegendreTri.h.

69 {
71 }
unsigned getIntegrationOrder() const
return current integration order.
MathLib::WeightedPoint getWeightedPoint(unsigned const igp) const

References getIntegrationOrder(), and getWeightedPoint().

Referenced by getWeightedPoint().

◆ getWeightedPoint() [3/3]

static MathLib::WeightedPoint NumLib::IntegrationGaussLegendreTri::getWeightedPoint ( unsigned const order,
unsigned const igp )
inlinestatic

get coordinates of a integration point

Parameters
orderthe number of integration points
igpthe sampling point id
Returns
weight

Definition at line 80 of file IntegrationGaussLegendreTri.h.

81 {
82 switch (order)
83 {
84 case 1:
85 return getWeightedPoint<MathLib::GaussLegendreTri<1>>(igp);
86 case 2:
87 return getWeightedPoint<MathLib::GaussLegendreTri<2>>(igp);
88 case 3:
89 return getWeightedPoint<MathLib::GaussLegendreTri<3>>(igp);
90 case 4:
91 return getWeightedPoint<MathLib::GaussLegendreTri<4>>(igp);
92 }
93 OGS_FATAL("Integration order {:d} not implemented for triangles.",
94 order);
95 }

References OGS_FATAL.

◆ setIntegrationOrder()

void NumLib::IntegrationGaussLegendreTri::setIntegrationOrder ( unsigned order)
inline

Change the integration order.

Definition at line 50 of file IntegrationGaussLegendreTri.h.

51 {
52 _order = order;
54 }
unsigned getNumberOfPoints() const
return the number of sampling points

References _n_sampl_pt, _order, and getNumberOfPoints().

Referenced by IntegrationGaussLegendreTri().

Member Data Documentation

◆ _n_sampl_pt

unsigned NumLib::IntegrationGaussLegendreTri::_n_sampl_pt {0}
private

Definition at line 128 of file IntegrationGaussLegendreTri.h.

128{0};

Referenced by getNumberOfPoints(), and setIntegrationOrder().

◆ _order

unsigned NumLib::IntegrationGaussLegendreTri::_order
private

Definition at line 127 of file IntegrationGaussLegendreTri.h.

Referenced by getIntegrationOrder(), and setIntegrationOrder().


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