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 Types

using WeightedPoint = MathLib::TemplateWeightedPoint< double, double, 2 >
 

Public Member Functions

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

Static Public Member Functions

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

Private Attributes

unsigned _order
 
unsigned _n_sampl_pt {0}
 

Member Typedef Documentation

◆ WeightedPoint

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 46 of file IntegrationGaussLegendreTri.h.

46  : _order(order)
47  {
48  this->setIntegrationOrder(order);
49  }
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 59 of file IntegrationGaussLegendreTri.h.

59 { 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 62 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 110 of file IntegrationGaussLegendreTri.h.

111  {
112  switch (order)
113  {
114  case 1:
116  case 2:
118  case 3:
120  case 4:
122  }
123  return 0;
124  }

◆ getWeightedPoint() [1/3]

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

Definition at line 99 of file IntegrationGaussLegendreTri.h.

100  {
101  return WeightedPoint(Method::X[igp], 0.5 * Method::W[igp]);
102  }
MathLib::TemplateWeightedPoint< double, double, 2 > WeightedPoint

◆ getWeightedPoint() [2/3]

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

get coordinates of a integration point

Parameters
igpThe integration point index
Returns
a weighted point

Definition at line 70 of file IntegrationGaussLegendreTri.h.

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

References getIntegrationOrder().

◆ getWeightedPoint() [3/3]

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

get coordinates of a integration point

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

Definition at line 82 of file IntegrationGaussLegendreTri.h.

83  {
84  switch (order)
85  {
86  case 1:
87  return getWeightedPoint<MathLib::GaussLegendreTri<1>>(igp);
88  case 2:
89  return getWeightedPoint<MathLib::GaussLegendreTri<2>>(igp);
90  case 3:
91  return getWeightedPoint<MathLib::GaussLegendreTri<3>>(igp);
92  case 4:
93  return getWeightedPoint<MathLib::GaussLegendreTri<4>>(igp);
94  }
95  return WeightedPoint(std::array<double, 2>(), 0);
96  }

◆ setIntegrationOrder()

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

Change the integration order.

Definition at line 52 of file IntegrationGaussLegendreTri.h.

53  {
54  _order = order;
56  }
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.

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: