Geophysical Inversion and Modelling Library v1.5.4
|
Public Member Functions | |
const R3Vector & | gauAbscissa (Index order) const |
const RVector & | gauWeights (Index order) const |
const R3Vector & | triGLAbscissa (Index order) const |
const RVector & | triGLWeights (Index order) const |
const R3Vector & | edgAbscissa (Index order) const |
const RVector & | edgWeights (Index order) const |
const R3Vector & | triAbscissa (Index order) const |
const RVector & | triWeights (Index order) const |
const R3Vector & | tetAbscissa (Index order) const |
const RVector & | tetWeights (Index order) const |
const R3Vector & | quaAbscissa (Index order) const |
const RVector & | quaWeights (Index order) const |
const R3Vector & | hexAbscissa (Index order) const |
const RVector & | hexWeights (Index order) const |
const R3Vector & | priAbscissa (Index order) const |
const RVector & | priWeights (Index order) const |
const R3Vector & | abscissa (const Shape &shape, uint order) const |
const RVector & | weights (const Shape &shape, uint order) const |
void | setTriUseGaussLegendre (bool use) |
bool | triUseGaussLegendre () const |
Protected Member Functions | |
void | initGau_ () |
void | initTriGL_ () |
void | initEdg_ () |
void | initTri_ () |
void | initTet_ () |
void | initQua_ () |
void | initHex_ () |
void | initPri_ () |
![]() | |
Singleton () | |
Singleton (const Singleton &) | |
Protected Attributes | |
bool | triUseGaussLegendre_ |
std::vector< R3Vector > | gauAbscissa_ |
std::vector< RVector > | gauWeights_ |
std::vector< R3Vector > | triGLAbscissa_ |
std::vector< RVector > | triGLWeights_ |
std::vector< R3Vector > | edgAbscissa_ |
std::vector< RVector > | edgWeights_ |
std::vector< R3Vector > | triAbscissa_ |
std::vector< RVector > | triWeights_ |
std::vector< R3Vector > | tetAbscissa_ |
std::vector< RVector > | tetWeights_ |
std::vector< R3Vector > | quaAbscissa_ |
std::vector< RVector > | quaWeights_ |
std::vector< R3Vector > | hexAbscissa_ |
std::vector< RVector > | hexWeights_ |
std::vector< R3Vector > | priAbscissa_ |
std::vector< RVector > | priWeights_ |
Friends | |
class | Singleton< IntegrationRules > |
Additional Inherited Members | |
![]() | |
static IntegrationRules * | pInstance () |
static IntegrationRules & | instance () |
Return Gauss-Legendre quadrature positions for a given shape of the MeshEntity upto order 10
References gauAbscissa(), GIMLI::Shape::rtti(), and triGLAbscissa().
const R3Vector & GIMLI::IntegrationRules::gauAbscissa | ( | Index | order | ) | const |
Return Gauss-Legendre quadrature point upto order <10.
Referenced by abscissa().
const RVector & GIMLI::IntegrationRules::gauWeights | ( | Index | order | ) | const |
Return Gauss-Legendre quadrature weights upto order <10.
Referenced by weights().
|
inline |
Set whether triangle integration use Gauss Legendre polynomials (up to order 9) or native triangle coordinates (up to order 5). Default is GL == true.
const R3Vector & GIMLI::IntegrationRules::triGLAbscissa | ( | Index | order | ) | const |
Generic quadrature positions for a triangle based on Gauss-Legendre quadrature H. T. RATHOD1*, K. V. NAGARAJA2, B. VENKATESUDU3 AND N. L. RAMESH4. Gauss Legendre quadrature over a triangle. J. Indian Inst. Sci., Sept.-Oct. 2004, 84, 183-188
References triGLAbscissa_.
Referenced by abscissa().
const RVector & GIMLI::IntegrationRules::triGLWeights | ( | Index | order | ) | const |
Generic quadrature weights for a triangle based on Gauss-Legendre quadrature H.T. RATHOD, K. V. NAGARAJA, B. VENKATESUDU AND N. L. RAMESH. Gauss Legendre quadrature over a triangle. J. Indian Inst. Sci., Sept.-Oct. 2004, 84, 183-188
References triGLWeights_.
Referenced by weights().
Return Gauss-Legendre quadrature weights for a given shape of the MeshEntity upto order 10
References gauWeights(), GIMLI::Shape::rtti(), and triGLWeights().
|
protected |
Quadrature positions for a triangle based on Gauss-Legendre quadrature
Referenced by triGLAbscissa().
|
protected |
Quadrature weight for a triangle based on Gauss-Legendre quadrature
Referenced by triGLWeights().