Class documentation of Concepts

Loading...
Searching...
No Matches

#include <neumannTrace.hh>

Inheritance diagram for hp2D::NTElement_BA< F >:
concepts::ElementWithCell< F > hp1D::IntegrableElm concepts::Element< F > concepts::IntegrationCell concepts::OutputOperator

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::Edgesupport () const
 
virtual concepts::Real3d vertex (uint i) const
 
virtual const concepts::Edge2dcell () 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::QuadratureRule1dintegration () 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::QuadRuleFactoryrule ()
 

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::QuadratureRule1dint_
 The integration rule.
 

Static Protected Attributes

static concepts::QuadRuleFactory rule_
 

Detailed Description

template<class F = Real>
class hp2D::NTElement_BA< F >

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.

Author
Kersten Schmidt, 2010, Robert Gruhlke 2012

Definition at line 39 of file neumannTrace.hh.

Member Typedef Documentation

◆ type

template<typename F >
typedef F concepts::ElementWithCell< F >::type
inherited

Definition at line 81 of file element.hh.

◆ UnderlyingElement

template<class F = Real>
typedef concepts::ElementAndFacette<hp2D::Element<F> > hp2D::NTElement_BA< F >::UnderlyingElement

Definition at line 42 of file neumannTrace.hh.

Member Enumeration Documentation

◆ intFormType

Integration form, which determines terms coming from integration over reference element

Definition at line 29 of file integral.hh.

Constructor & Destructor Documentation

◆ NTElement_BA()

template<class F = Real>
hp2D::NTElement_BA< F >::NTElement_BA ( const concepts::Edge2d cell)

Constructor

Creates the Neumann trace from the element on one side.

Member Function Documentation

◆ addElement()

template<class F = Real>
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.

◆ appendT()

template<class F = Real>
void hp2D::NTElement_BA< F >::appendT ( concepts::TColumn< F > *  T)
inline

Appends the T columns to the T matrix.

Definition at line 123 of file neumannTrace.hh.

◆ cell()

template<class F = Real>
virtual const concepts::Edge2d & hp2D::NTElement_BA< F >::cell ( ) const
inlinevirtual

Returns the cell on which the element is built.

Implements concepts::ElementWithCell< F >.

Definition at line 76 of file neumannTrace.hh.

◆ chi()

const concepts::Real3d hp1D::IntegrableElm::chi ( const Real  x) const
inlineinherited

Computes the element map. The reference element is [0,1].

Definition at line 40 of file element.hh.

◆ elemMap() [1/3]

template<typename F >
Real3d concepts::ElementWithCell< F >::elemMap ( const Real  coord_local) const
inlineinherited

Definition at line 86 of file element.hh.

◆ elemMap() [2/3]

template<typename F >
Real3d concepts::ElementWithCell< F >::elemMap ( const Real2d coord_local) const
inlineinherited

Definition at line 90 of file element.hh.

◆ elemMap() [3/3]

template<typename F >
Real3d concepts::ElementWithCell< F >::elemMap ( const Real3d coord_local) const
inlineinherited

Definition at line 94 of file element.hh.

◆ graphics()

◆ info()

template<class F = Real>
virtual std::ostream & hp2D::NTElement_BA< F >::info ( std::ostream &  os) const
protectedvirtual

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

◆ integration()

const concepts::QuadratureRule1d * hp1D::IntegrableElm::integration ( ) const
inlineinherited

Returns the integration rule.

Definition at line 51 of file element.hh.

◆ jacobianDeterminant()

Real hp1D::IntegrableElm::jacobianDeterminant ( const Real  x) const
inlineinherited

Computes the determinant of the Jacobian.

Definition at line 45 of file element.hh.

◆ mapRho()

template<class F = Real>
const concepts::Z2 hp2D::NTElement_BA< F >::mapRho ( ) const
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.

◆ p()

template<class F = Real>
ushort hp2D::NTElement_BA< F >::p ( ) const
inline

Definition at line 127 of file neumannTrace.hh.

◆ quadraturePoint()

virtual bool hp1D::IntegrableElm::quadraturePoint ( uint  i,
intPoint p,
intFormType  form = ZERO,
bool  localCoord = false 
) const
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.

Parameters
inumber of quadrature point
intPointdata given back
formIntegration form
localCoordIf true, local coordinates are returned. Else physical coordinates.

Implements concepts::IntegrationCell.

◆ recomputeShapefunctions()

template<class F = Real>
void hp2D::NTElement_BA< F >::recomputeShapefunctions ( )

Recompute shape functions, e.g. for other abscissas redefined through IntegrableElm::rule().set(...)

◆ rule()

static concepts::QuadRuleFactory & hp1D::IntegrableElm::rule ( )
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.

◆ shpfct()

template<class F = Real>
const concepts::ShapeFunction1D< Real > * hp2D::NTElement_BA< F >::shpfct ( ) const
inline

Returns the shape functions.

Definition at line 86 of file neumannTrace.hh.

◆ shpfctD()

template<class F = Real>
const concepts::ShapeFunction1D< Real > * hp2D::NTElement_BA< F >::shpfctD ( ) const
inline

Returns the derivatives of the shape functions

TODO: implementation.

Definition at line 92 of file neumannTrace.hh.

◆ support()

template<class F = Real>
virtual const concepts::Edge & hp2D::NTElement_BA< F >::support ( ) const
inlinevirtual

Definition at line 73 of file neumannTrace.hh.

◆ T()

template<class F = Real>
virtual const concepts::TMatrix< F > & hp2D::NTElement_BA< F >::T ( ) const
inlinevirtual

Returns the T matrix of the element.

Implements concepts::ElementWithCell< F >.

Definition at line 78 of file neumannTrace.hh.

◆ tag()

template<class F >
uint & concepts::Element< F >::tag ( )
inlineinherited

Returns the tag.

Definition at line 66 of file element.hh.

◆ uelm()

template<class F = Real>
const concepts::Sequence< UnderlyingElement > & hp2D::NTElement_BA< F >::uelm ( ) const
inline

Returns the Underlying Element(s)

Definition at line 101 of file neumannTrace.hh.

Member Data Documentation

◆ int_

std::unique_ptr<concepts::QuadratureRule1d> hp1D::IntegrableElm::int_
protectedinherited

The integration rule.

Definition at line 72 of file element.hh.

◆ rule_

concepts::QuadRuleFactory hp1D::IntegrableElm::rule_
staticprotectedinherited

Definition at line 74 of file element.hh.

◆ T_

template<class F = Real>
concepts::TMatrix<Real> hp2D::NTElement_BA< F >::T_
protected

T matrix of the element.

Definition at line 139 of file neumannTrace.hh.


The documentation for this class was generated from the following file: