OGS
NumLib::IntegrationGaussLegendreRegular< N_DIM > Class Template Reference

Detailed Description

template<unsigned N_DIM>
class NumLib::IntegrationGaussLegendreRegular< N_DIM >

Gauss-Legendre quadrature rule for regular shape elements: line, quad and hex.

Template Parameters
N_DIMSpatial dimension

Definition at line 29 of file IntegrationGaussLegendreRegular.h.

#include <IntegrationGaussLegendreRegular.h>

Public Member Functions

 IntegrationGaussLegendreRegular (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
 
std::array< unsigned, 1 > getPositionIndices (unsigned, unsigned igp)
 
std::array< unsigned, 2 > getPositionIndices (unsigned order, unsigned igp)
 
std::array< unsigned, 3 > getPositionIndices (unsigned order, unsigned igp)
 

Static Public Member Functions

static std::array< unsigned, N_DIM > getPositionIndices (unsigned order, unsigned igp)
 
static MathLib::WeightedPoint getWeightedPoint (unsigned const order, unsigned const igp)
 

Static Private Member Functions

template<typename Method >
static MathLib::WeightedPoint getWeightedPoint (std::array< unsigned, N_DIM > const &pos)
 

Private Attributes

unsigned _order
 
unsigned _n_sampl_pt {0}
 

Constructor & Destructor Documentation

◆ IntegrationGaussLegendreRegular()

template<unsigned N_DIM>
NumLib::IntegrationGaussLegendreRegular< N_DIM >::IntegrationGaussLegendreRegular ( unsigned order = 2)
inlineexplicit

Create IntegrationGaussLegendreRegular of the given Gauss-Legendre integration order.

Parameters
orderintegration order (default 2)

Definition at line 36 of file IntegrationGaussLegendreRegular.h.

36 : _order(order)
37 {
38 this->setIntegrationOrder(order);
39 }
void setIntegrationOrder(unsigned order)
Change the integration order.

References NumLib::IntegrationGaussLegendreRegular< N_DIM >::setIntegrationOrder().

Member Function Documentation

◆ getIntegrationOrder()

template<unsigned N_DIM>
unsigned NumLib::IntegrationGaussLegendreRegular< N_DIM >::getIntegrationOrder ( ) const
inline

◆ getNumberOfPoints()

template<unsigned N_DIM>
unsigned NumLib::IntegrationGaussLegendreRegular< N_DIM >::getNumberOfPoints ( ) const
inline

◆ getPositionIndices() [1/4]

std::array< unsigned, 2 > NumLib::IntegrationGaussLegendreRegular< 2 >::getPositionIndices ( unsigned order,
unsigned igp )
inline

Definition at line 29 of file IntegrationGaussLegendreRegular-impl.h.

31{
32 assert(igp < order * order);
33 std::array<unsigned, 2> result;
34 result[0] = igp / order;
35 result[1] = igp % order;
36 return result;
37}

◆ getPositionIndices() [2/4]

std::array< unsigned, 3 > NumLib::IntegrationGaussLegendreRegular< 3 >::getPositionIndices ( unsigned order,
unsigned igp )
inline

Definition at line 41 of file IntegrationGaussLegendreRegular-impl.h.

43{
44 assert(igp < order * order * order);
45 unsigned const gp_r = igp / (order * order);
46 unsigned const gp_s = igp % (order * order);
47 std::array<unsigned, 3> result;
48 result[0] = gp_r;
49 result[1] = gp_s / order;
50 result[2] = gp_s % order;
51 return result;
52}

◆ getPositionIndices() [3/4]

template<unsigned N_DIM>
static std::array< unsigned, N_DIM > NumLib::IntegrationGaussLegendreRegular< N_DIM >::getPositionIndices ( unsigned order,
unsigned igp )
static

Get position indexes in r-s-t axis.

Parameters
orderThe number of integration points
igpThe integration point index
Returns
a tuple of position indexes

◆ getPositionIndices() [4/4]

std::array< unsigned, 1 > NumLib::IntegrationGaussLegendreRegular< 1 >::getPositionIndices ( unsigned ,
unsigned igp )
inline

Definition at line 19 of file IntegrationGaussLegendreRegular-impl.h.

21{
22 std::array<unsigned, 1> result;
23 result[0] = igp;
24 return result;
25}

◆ getWeightedPoint() [1/3]

template<unsigned N_DIM>
template<typename Method >
MathLib::WeightedPoint NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint ( std::array< unsigned, N_DIM > const & pos)
inlinestaticprivate

Computes weighted point using given integration method.

Template Parameters
MethodIntegration method to use.
Parameters
posPoint indices computed by getPositionIndices.

Definition at line 79 of file IntegrationGaussLegendreRegular-impl.h.

81{
82 std::array<double, N_DIM> coords;
83 double weight = 1;
84 for (unsigned d = 0; d < N_DIM; d++)
85 {
86 coords[d] = Method::X[pos[d]];
87 weight *= Method::W[pos[d]];
88 }
89
90 return MathLib::WeightedPoint(coords, weight);
91}
constexpr ranges::views::view_closure coords
Definition Mesh.h:232

◆ getWeightedPoint() [2/3]

template<unsigned N_DIM>
MathLib::WeightedPoint NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint ( unsigned const igp) const
inline

Get coordinates of the integration point.

Parameters
igpThe integration point index
Returns
a weighted point

Definition at line 58 of file IntegrationGaussLegendreRegular.h.

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

References NumLib::IntegrationGaussLegendreRegular< N_DIM >::getIntegrationOrder(), and NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint().

Referenced by NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint().

◆ getWeightedPoint() [3/3]

template<unsigned N_DIM>
MathLib::WeightedPoint NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint ( unsigned const order,
unsigned const igp )
inlinestatic

Get coordinates of the integration point.

Parameters
orderThe number of integration points
igpThe integration point index
Returns
a weighted point

Definition at line 56 of file IntegrationGaussLegendreRegular-impl.h.

58{
59 assert(igp < std::pow(order, N_DIM));
60 std::array<unsigned, N_DIM> const pos = getPositionIndices(order, igp);
61
62 switch (order)
63 {
64 case 1:
65 return getWeightedPoint<MathLib::GaussLegendre<1>>(pos);
66 case 2:
67 return getWeightedPoint<MathLib::GaussLegendre<2>>(pos);
68 case 3:
69 return getWeightedPoint<MathLib::GaussLegendre<3>>(pos);
70 case 4:
71 return getWeightedPoint<MathLib::GaussLegendre<4>>(pos);
72 }
73 OGS_FATAL("Integration order {:d} not implemented.", order);
74}
#define OGS_FATAL(...)
Definition Error.h:26
static std::array< unsigned, N_DIM > getPositionIndices(unsigned order, unsigned igp)

References OGS_FATAL.

◆ setIntegrationOrder()

template<unsigned N_DIM>
void NumLib::IntegrationGaussLegendreRegular< N_DIM >::setIntegrationOrder ( unsigned order)
inline

Change the integration order.

Definition at line 42 of file IntegrationGaussLegendreRegular.h.

43 {
44 this->_n_sampl_pt = static_cast<unsigned>(std::pow(order, N_DIM));
45 this->_order = order;
46 }

References NumLib::IntegrationGaussLegendreRegular< N_DIM >::_n_sampl_pt, and NumLib::IntegrationGaussLegendreRegular< N_DIM >::_order.

Referenced by NumLib::IntegrationGaussLegendreRegular< N_DIM >::IntegrationGaussLegendreRegular().

Member Data Documentation

◆ _n_sampl_pt

◆ _order


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