#include <quad.hh>
Public Types | |
typedef F | type |
enum | intFormType { ZERO , ONE , TWO , THREE } |
Public Member Functions | |
BaseQuad (concepts::QuadNd &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1) | |
virtual const concepts::Quad & | support () const |
virtual concepts::Real3d | vertex (uint i) const |
virtual const concepts::QuadNd & | cell () const |
Returns the cell on which the element is built. | |
virtual const concepts::ElementGraphics< F > * | graphics () const =0 |
Returns element graphics class. | |
virtual const concepts::TMatrix< F > & | T () const |
Returns the T matrix of the element. | |
void | appendT (concepts::TColumn< F > *T) |
Appends the T columns to the T matrix. | |
Real3d | elemMap (const Real coord_local) const |
Real3d | elemMap (const Real2d &coord_local) const |
Real3d | elemMap (const Real3d &coord_local) const |
uint & | tag () |
Returns the tag. | |
concepts::Real2d | chi (const Real x, const Real y) const |
concepts::MapReal2d | jacobian (const Real x, const Real y) const |
Computes the Jacobian. | |
concepts::MapReal2d | jacobianInverse (const Real x, const Real y) const |
Computes the inverse of the Jacobian. | |
concepts::MapReal2d | inverseLaplace (const Real x, const Real y) const |
Real | gramDeterminantRoot (const Real x, const Real y) const |
void | setStrategy (const concepts::Quad2dSubdivision *strategy=0) |
const concepts::Quad2dSubdivision * | getStrategy () const |
virtual bool | quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const |
const concepts::QuadratureRule2d *const | intRule () const |
Static Public Member Functions | |
static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & | factory () |
static concepts::QuadRuleFactoryTensor2d * | factory_rp () |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. | |
Protected Attributes | |
concepts::TMatrix< F > | T_ |
T matrix of the element. | |
std::unique_ptr< concepts::QuadratureRule2d > | intRule_ |
The integration rules. | |
concepts::QuadNd & | cell_ |
The cell. | |
Static Protected Attributes | |
static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > | factory_ |
A base of a 2D quad FEM element for different basis functions
|
inherited |
Definition at line 81 of file element.hh.
|
inherited |
Integration form, which determines terms coming from integration over reference element
Definition at line 29 of file integral.hh.
hp2D::BaseQuad< F >::BaseQuad | ( | concepts::QuadNd & | cell, |
const ushort * | p, | ||
concepts::TColumn< F > * | T0, | ||
concepts::TColumn< F > * | T1 | ||
) |
Constructor
cell | Cell on which the element is defined |
p | Polynomial degree (might be anisotropic and is only) taking to initialize the quadrature rule |
T0 | Part of the T matrix |
T1 | Part of the T matrix |
|
inlineinherited |
Appends the T columns to the T matrix.
Definition at line 43 of file element.hh.
|
inlinevirtual |
Returns the cell on which the element is built.
Implements hp2D::Element< F >.
|
inlineinherited |
|
inlineinherited |
Definition at line 86 of file element.hh.
|
inlineinherited |
Definition at line 90 of file element.hh.
|
inlineinherited |
Definition at line 94 of file element.hh.
|
inlinestaticinherited |
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.
|
inlinestaticinherited |
|
inlineinherited |
|
inlineinherited |
|
pure virtual |
Returns element graphics class.
Reimplemented from concepts::Element< F >.
Implemented in hp2D::Quad< F >, hp2D::Quad< H >, hp2D::Quad< Real >, and hp2Dedge::Quad< F >.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from hp2D::Element< F >.
Reimplemented in hp2D::Quad< F >, hp2D::Quad< H >, hp2D::Quad< Real >, hp2Dedge::Quad< F >, and hp3D::NeumannTraceElement3d< F >.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
virtualinherited |
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.
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.
|
inlineinherited |
Sets the subdivision strategy of the underlying cell of this element. It calls Quad2d::setStrategy.
strategy | Pointer to an instance of a subdivision strategy. |
StrategyChange | if the change is not allowed (the change is not allowed if there are children present) |
|
inlinevirtual |
Returns the topological support of the element. Possible supports for an element are quadrilaterals and triangles.
Implements hp2D::Element< F >.
|
inlinevirtualinherited |
Returns the T matrix of the element.
Implements concepts::ElementWithCell< F >.
Reimplemented in hp3D::NeumannTraceElement3d< F >.
Definition at line 40 of file element.hh.
|
inlineinherited |
Returns the tag.
Definition at line 66 of file element.hh.
|
inlinevirtual |
Returns the coordinates of the ith vertex of this element.
i | Index of the vertex |
Implements hp2D::Element< F >.
|
protectedinherited |
|
staticprotectedinherited |
|
protectedinherited |
|
protectedinherited |
T matrix of the element.
Definition at line 57 of file element.hh.