hp2D::IntegrableQuad Class Reference
Class holding the quadrature rule and the cell of a quadrilateral element. More...
#include <quad.hh>
Public Types | |
enum | intFormType { ZERO, ONE, TWO, THREE } |
Integration form, which determines terms coming from integration over reference element. More... | |
Public Member Functions | |
concepts::Real2d | chi (const Real x, const Real y) const |
Computes the element map, the map from the reference element to the real geometry of the quad. More... | |
const concepts::Quad2dSubdivision * | getStrategy () const |
Returns the subdivision strategy of the underlying cell of this element. More... | |
Real | gramDeterminantRoot (const Real x, const Real y) const |
Computes the sqaure root of the Gram determinant multiplied by 0.25 to incorporate the transformation factor for the integration weights. More... | |
IntegrableQuad (concepts::QuadNd &cell) | |
const concepts::QuadratureRule2d *const | intRule () const |
concepts::MapReal2d | inverseLaplace (const Real x, const Real y) const |
concepts::MapReal2d | jacobian (const Real x, const Real y) const |
Computes the Jacobian. More... | |
concepts::MapReal2d | jacobianInverse (const Real x, const Real y) const |
Computes the inverse of the Jacobian. More... | |
virtual bool | quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const |
Delivers a quadrature point. More... | |
void | setStrategy (const concepts::Quad2dSubdivision *strategy=0) throw (concepts::StrategyChange) |
Sets the subdivision strategy of the underlying cell of this element. More... | |
Static Public Member Functions | |
static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & | factory () |
Access to the quadrature rule, which is valid for all elements of this type (hp2D::IntegrableQuad). More... | |
static concepts::QuadRuleFactoryTensor2d * | factory_rp () |
Access to the quadrature rule ( as factory() ) through a raw pointer. More... | |
Protected Attributes | |
concepts::QuadNd & | cell_ |
The cell. More... | |
std::unique_ptr< concepts::QuadratureRule2d > | intRule_ |
The integration rules. More... | |
Static Protected Attributes | |
static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > | factory_ |
Detailed Description
Class holding the quadrature rule and the cell of a quadrilateral element.
Member Enumeration Documentation
◆ intFormType
|
inherited |
Integration form, which determines terms coming from integration over reference element.
Enumerator | |
---|---|
ZERO | |
ONE | |
TWO | |
THREE |
Definition at line 29 of file integral.hh.
Constructor & Destructor Documentation
◆ IntegrableQuad()
hp2D::IntegrableQuad::IntegrableQuad | ( | concepts::QuadNd & | cell | ) |
Member Function Documentation
◆ chi()
|
inline |
◆ factory()
|
inlinestatic |
Access to the quadrature rule, which is valid for all elements of this type (hp2D::IntegrableQuad).
Change of the quadrature rule is put into practice for newly created elements and for already created elements by precomputing the integration points and shape functions on them.
◆ factory_rp()
|
inlinestatic |
◆ getStrategy()
|
inline |
◆ gramDeterminantRoot()
|
inline |
◆ intRule()
|
inline |
◆ inverseLaplace()
|
inline |
◆ jacobian()
|
inline |
◆ jacobianInverse()
|
inline |
◆ quadraturePoint()
|
virtual |
Delivers a quadrature point.
Quadrature point consists of coordinates (for evaluation of formulas) and intermediate data, consisting of the weight and term coming from mapping.
Returns false
, if the number of quadrature points is overstepped.
- Parameters
-
i number of quadrature point intPoint data given back form Integration form localCoord If true, local coordinates are returned. Else physical coordinates.
Implements concepts::IntegrationCell.
◆ setStrategy()
|
inline |
Sets the subdivision strategy of the underlying cell of this element.
It calls Quad2d::setStrategy.
- Parameters
-
strategy Pointer to an instance of a subdivision strategy.
- Exceptions
-
StrategyChange if the change is not allowed (the change is not allowed if there are children present)
Member Data Documentation
◆ cell_
|
protected |
◆ factory_
|
staticprotected |
◆ intRule_
|
protected |
The documentation for this class was generated from the following file:
- hp2D/quad.hh