Loading [MathJax]/extensions/tex2jax.js
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 28 of file IntegrationGaussLegendreRegular.h.

#include <IntegrationGaussLegendreRegular.h>

Public Types

using WeightedPoint = typename MathLib::TemplateWeightedPoint< double, double, N_DIM >
 

Public Member Functions

 IntegrationGaussLegendreRegular (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
 
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)
 
template<typename Method >
MathLib::TemplateWeightedPoint< double, double, N_DIM > getWeightedPoint (std::array< unsigned, N_DIM > const &pos)
 

Static Public Member Functions

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

Static Private Member Functions

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

Private Attributes

unsigned _order
 
unsigned _n_sampl_pt {0}
 

Member Typedef Documentation

◆ WeightedPoint

template<unsigned N_DIM>
using NumLib::IntegrationGaussLegendreRegular< N_DIM >::WeightedPoint = typename MathLib::TemplateWeightedPoint<double, double, N_DIM>

Definition at line 31 of file IntegrationGaussLegendreRegular.h.

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 38 of file IntegrationGaussLegendreRegular.h.

38  : _order(order)
39  {
40  this->setIntegrationOrder(order);
41  }
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

return current integration order.

Definition at line 51 of file IntegrationGaussLegendreRegular.h.

51 { return _order; }

References NumLib::IntegrationGaussLegendreRegular< N_DIM >::_order.

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

◆ getNumberOfPoints()

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

return the number of sampling points

Definition at line 54 of file IntegrationGaussLegendreRegular.h.

References NumLib::IntegrationGaussLegendreRegular< N_DIM >::_n_sampl_pt.

◆ 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/4]

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

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

83 {
84  std::array<double, N_DIM> coords;
85  double weight = 1;
86  for (unsigned d = 0; d < N_DIM; d++)
87  {
88  coords[d] = Method::X[pos[d]];
89  weight *= Method::W[pos[d]];
90  }
91 
93  weight);
94 }

◆ getWeightedPoint() [2/4]

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

Computes weighted point using given integration method.

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

◆ getWeightedPoint() [3/4]

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

Get coordinates of the integration point.

Parameters
igpThe integration point index
Returns
a weighted point

Definition at line 60 of file IntegrationGaussLegendreRegular.h.

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

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

◆ getWeightedPoint() [4/4]

template<unsigned N_DIM>
MathLib::TemplateWeightedPoint< double, double, N_DIM > NumLib::IntegrationGaussLegendreRegular< N_DIM >::getWeightedPoint ( unsigned  order,
unsigned  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 
75  std::array<double, N_DIM>(), 0);
76 }
static std::array< unsigned, N_DIM > getPositionIndices(unsigned order, unsigned igp)

◆ setIntegrationOrder()

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

Change the integration order.

Definition at line 44 of file IntegrationGaussLegendreRegular.h.

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

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: