#include <neumannTrace.hh>
Classes | |
class | ShapeFunction |
Public Types | |
typedef concepts::ElementAndFacette< hp2D::Element< F > > | UnderlyingElement |
typedef F | type |
enum | intFormType { ZERO , ONE , TWO , THREE } |
Public Member Functions | |
NTElement_BA (const concepts::Edge2d &cell) | |
virtual const concepts::Edge & | support () const |
virtual concepts::Real3d | vertex (uint i) const |
virtual const concepts::Edge2d & | cell () const |
Returns the cell on which the element is built. | |
virtual const concepts::TMatrix< F > & | T () const |
Returns the T matrix of the element. | |
void | recomputeShapefunctions () |
const concepts::ShapeFunction1D< Real > * | shpfct () const |
Returns the shape functions. | |
const concepts::ShapeFunction1D< Real > * | shpfctD () const |
const concepts::Sequence< UnderlyingElement > & | uelm () const |
const concepts::Z2 | mapRho () const |
void | addElement (const hp2D::Quad< Real > &quad, uint k, Real weight=1.0) |
void | appendT (concepts::TColumn< F > *T) |
Appends the T columns to the T matrix. | |
ushort | p () const |
Real3d | elemMap (const Real coord_local) const |
Real3d | elemMap (const Real2d &coord_local) const |
Real3d | elemMap (const Real3d &coord_local) const |
virtual const ElementGraphics< F > * | graphics () 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 | |
concepts::TMatrix< Real > | T_ |
T matrix of the element. | |
std::unique_ptr< concepts::QuadratureRule1d > | int_ |
The integration rule. | |
Static Protected Attributes | |
static concepts::QuadRuleFactory | rule_ |
Element on an edge representing the normal derivative of neighbouring elements, especially their mean or jump.
The number of shape functions are the sum of the number of shape functions of the two neighbouring 2D elements.
Definition at line 39 of file neumannTrace.hh.
|
inherited |
Definition at line 81 of file element.hh.
typedef concepts::ElementAndFacette<hp2D::Element<F> > hp2D::NTElement_BA< F >::UnderlyingElement |
Definition at line 42 of file neumannTrace.hh.
|
inherited |
Integration form, which determines terms coming from integration over reference element
Definition at line 29 of file integral.hh.
hp2D::NTElement_BA< F >::NTElement_BA | ( | const concepts::Edge2d & | cell | ) |
Constructor
Creates the Neumann trace from the element on one side.
void hp2D::NTElement_BA< 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.
|
inline |
Appends the T columns to the T matrix.
Definition at line 123 of file neumannTrace.hh.
|
inlinevirtual |
Returns the cell on which the element is built.
Implements concepts::ElementWithCell< F >.
Definition at line 76 of file neumannTrace.hh.
|
inlineinherited |
Computes the element map. The reference element is [0,1].
Definition at line 40 of file element.hh.
|
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.
|
inlinevirtualinherited |
Reimplemented in hp1D::BaseElement< F >, hp1D::BaseElement< Real >, hp2D::Quad< F >, hp2D::Quad< H >, hp2D::Quad< Real >, hp2Dedge::Quad< F >, linDG2D::Triangle, linDG3D::FvdgP0TetElem, linDG3D::FvdgP1TetElem, hp2D::BaseQuad< F >, linDG3D::FvdgElement, hp2D::BaseQuad< H >, and hp2D::BaseQuad< Real >.
Definition at line 63 of file element.hh.
|
protectedvirtual |
Returns information in an output stream.
Reimplemented from concepts::OutputOperator.
|
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.
|
inline |
Returns the orientation of the MappingEdge2d corresponding to the quad. Analogue to understand like topological orientations. I.e. this method just makes sense after adding a quad TODO: For 2 Quads make sequence<concepts::Z2> datatype returns : 0 if mapping is counter-clockwise 1 if mapping is clock-wise
Definition at line 113 of file neumannTrace.hh.
|
inline |
Definition at line 127 of file neumannTrace.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.
void hp2D::NTElement_BA< F >::recomputeShapefunctions | ( | ) |
Recompute shape functions, e.g. for other abscissas redefined through IntegrableElm::rule().set(...)
|
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.
|
inline |
Returns the shape functions.
Definition at line 86 of file neumannTrace.hh.
|
inline |
Returns the derivatives of the shape functions
TODO: implementation.
Definition at line 92 of file neumannTrace.hh.
|
inlinevirtual |
Definition at line 73 of file neumannTrace.hh.
|
inlinevirtual |
Returns the T matrix of the element.
Implements concepts::ElementWithCell< F >.
Definition at line 78 of file neumannTrace.hh.
|
inlineinherited |
Returns the tag.
Definition at line 66 of file element.hh.
|
inline |
Returns the Underlying Element(s)
Definition at line 101 of file neumannTrace.hh.
|
protectedinherited |
The integration rule.
Definition at line 72 of file element.hh.
|
staticprotectedinherited |
Definition at line 74 of file element.hh.
|
protected |
T matrix of the element.
Definition at line 139 of file neumannTrace.hh.