#include <neumannTraceElement.hh>
Classes | |
struct | mapPart |
class | NTShapeFunction |
Public Types | |
typedef concepts::ElementAndFacette< hp2D::Element< F > > | UnderlyingElement |
typedef F | FieldT |
typedef F | type |
enum | intFormType { ZERO , ONE , TWO , THREE } |
Public Member Functions | |
NeumannTraceElement (const concepts::EdgeNd &cell, uint p, const mapPart *coarseCell=0) | |
virtual const concepts::TMatrix< F > & | T () const |
Returns the T matrix of the element. | |
void | computeShpfctVals (const concepts::Array< Real > &xP, concepts::Array< Real > &val) const |
virtual void | recomputeShapefunctions () |
void | buildT () |
virtual const concepts::ShapeFunction1D< Real > * | shpfct () const |
Returns the shape functions. | |
const concepts::Sequence< UnderlyingElement > & | uelm () const |
void | addElement (const hp2D::Quad< Real > &quad, uint k, Real weight=1.0) |
void | addIrrElement (const hp2D::Quad< Real > &coarseQuad, uint k_coarse, Real weight) |
virtual const concepts::Edge & | support () const |
virtual concepts::Real3d | vertex (uint i) const |
virtual const concepts::EdgeNd & | cell () const |
Returns the cell on which the element is built. | |
ushort | p () const |
virtual const concepts::ElementGraphics< Real > * | graphics () const |
Returns element graphics class. | |
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. | |
const concepts::Real3d | chi (const Real x) const |
Real | jacobianDeterminant (const Real x) const |
Computes the determinant of the Jacobian. | |
const concepts::QuadratureRule1d * | integration () const |
Returns the integration rule. | |
virtual bool | quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const |
Static Public Member Functions | |
static concepts::QuadRuleFactory & | rule () |
Protected Member Functions | |
virtual std::ostream & | info (std::ostream &os) const |
Returns information in an output stream. | |
Protected Attributes | |
ushort | p_ |
T matrix of the element. | |
const concepts::EdgeNd & | cell_ |
The cell. | |
std::unique_ptr< concepts::QuadratureRule1d > | int_ |
The integration rule. | |
Static Protected Attributes | |
static std::unique_ptr< LineGraphics > | graphics_ |
static concepts::QuadRuleFactory | rule_ |
Element on an edge representing the normal derivatives of neighbouring elements, especially their mean or jump.
Also a PLUS and MINUS setting is involded, that is a FIRST trace kind, with respect to the K+ and K- introduced elements via the setted Normalvector.
The number of shape functions are the sum of the number of shape functions of the two neighbouring 2D elements.
Definition at line 42 of file neumannTraceElement.hh.
|
inherited |
Definition at line 82 of file element.hh.
|
inherited |
Definition at line 81 of file element.hh.
typedef concepts::ElementAndFacette<hp2D::Element<F> > hp2D::NeumannTraceElement< F >::UnderlyingElement |
Definition at line 44 of file neumannTraceElement.hh.
|
inherited |
Integration form, which determines terms coming from integration over reference element
Definition at line 29 of file integral.hh.
hp2D::NeumannTraceElement< F >::NeumannTraceElement | ( | const concepts::EdgeNd & | cell, |
uint | p, | ||
const mapPart * | coarseCell = 0 |
||
) |
Constructor
For irregular meshes the underlying elements of a neumanntraceelement is a fine Quad and a neighboured coarse Quad. The coarseCell describes the edgemapping of the coarseedge of which the cell is a subset (part). The coarseCellpart also marks if a neumanntraceElement is build upon regular or irregular mesh
cell | geometric cell of this element |
p | highest Polynomial degree of underlying Elements |
coarseCell | for irregular meshes where cell's underlying |
void hp2D::NeumannTraceElement< F >::addElement | ( | const hp2D::Quad< Real > & | quad, |
uint | k, | ||
Real | weight = 1.0 |
||
) |
Adds the contribution to the Neumann trace from the element on one of the (at most two) sides. Basis functions between quads do not have to be continous. In this case the highest polynomial degree along edge is set, i.e. for integration rule
quad | Quad with k-th edge is NeumannTraceElement's topological cell. |
void hp2D::NeumannTraceElement< F >::addIrrElement | ( | const hp2D::Quad< Real > & | coarseQuad, |
uint | k_coarse, | ||
Real | weight | ||
) |
Basically the same as addElement, but with the difference in the meaning. In this case the general ancestor edge of this finer Edge is the k-th edge of the coarseQuad.
coarseQuad | Quad with k-th edge is general ancestor of NeumannTraceElement's topological cell. |
void hp2D::NeumannTraceElement< F >::buildT | ( | ) |
After adding the Elements this routine builds the TMatrix and polynomial informations on the element.
|
inlinevirtualinherited |
Returns the cell on which the element is built.
Implements concepts::ElementWithCell< F >.
Definition at line 99 of file element.hh.
|
inlineinherited |
Computes the element map. The reference element is [0,1].
Definition at line 40 of file element.hh.
void hp2D::NeumannTraceElement< F >::computeShpfctVals | ( | const concepts::Array< Real > & | xP, |
concepts::Array< Real > & | val | ||
) | const |
Evaluation method to get basisfunctions evaluated on requested points this should be used for functions, i.e. hp2D::Value
xP | requested Points to evaluate on |
val | evaluated shpfct, method resizes array. |
|
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.
|
virtualinherited |
Returns element graphics class.
Reimplemented from concepts::Element< F >.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from hp1D::BaseElement< F >.
|
inlineinherited |
Returns the integration rule.
Definition at line 51 of file element.hh.
|
inlineinherited |
Computes the determinant of the Jacobian.
Definition at line 45 of file element.hh.
|
inlineinherited |
Definition at line 103 of file element.hh.
|
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.
|
virtual |
Recompute shape functions, e.g. for other abscissas redefined through IntegrableElm::rule().set(...)
Implements hp1D::BaseElement< F >.
|
inlinestaticinherited |
Access to the quadrature rule, which is valid for all elements of this type (hp1D::IntegrableElm).
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.
Definition at line 62 of file element.hh.
|
virtual |
Returns the shape functions.
Implements hp1D::BaseElement< F >.
|
inlinevirtualinherited |
Definition at line 91 of file element.hh.
|
inlinevirtual |
Returns the T matrix of the element.
Implements concepts::ElementWithCell< F >.
Definition at line 121 of file neumannTraceElement.hh.
|
inlineinherited |
Returns the tag.
Definition at line 66 of file element.hh.
|
inline |
Returns the Underlying Element(s)
Definition at line 151 of file neumannTraceElement.hh.
|
inlinevirtualinherited |
Definition at line 95 of file element.hh.
|
protectedinherited |
The cell.
Definition at line 70 of file element.hh.
|
staticprotectedinherited |
Definition at line 126 of file element.hh.
|
protectedinherited |
The integration rule.
Definition at line 72 of file element.hh.
|
protectedinherited |
T matrix of the element.
Definition at line 124 of file element.hh.
|
staticprotectedinherited |
Definition at line 74 of file element.hh.