Class documentation of Concepts

Loading...
Searching...
No Matches
hp2Dedge::Edge< F > Class Template Reference

#include <quad.hh>

Inheritance diagram for hp2Dedge::Edge< F >:
concepts::Element< F > concepts::OutputOperator

Public Types

typedef F type
 

Public Member Functions

 Edge (const Quad< F > &elm, const ushort k)
 
concepts::Real2d vertex (uint i) const
 
virtual const concepts::Edgesupport () const
 
const ushort edge () const
 Returns number of the edge.
 
const Quad< F > & elm () const
 Returns element.
 
const ushort direction () const
 Returns direction of edge on reference quad [0,1]^2, 0 - x, 1 - y.
 
const Real sign () const
 Returns sign of outer normal vector, e.g. left edge -1, right edge +1.
 
const hp2D::KarniadakisDeriv2shpfct () const
 
const concepts::QuadratureRule1dintegration () const
 Returns the integration rule.
 
concepts::Real2d localCoords (const Real t) const
 coordinate of point on the edge inside reference element [0,1]^2
 
concepts::Real2d chi (const Real t) const
 
concepts::MapReal2d jacobian (const Real t) const
 
concepts::MapReal2d jacobianInverse (const Real t) const
 Computes the inverse of the Jacobian.
 
Real jacobianDeterminant (const Real t) const
 Computes the determinant of the Jacobian.
 
Real diffElement (const Real t) const
 Computes the differential element for integration over [-1,1].
 
virtual const concepts::TMatrixBase< F > & T () const
 T-Matrix of the appropiate Quad, not used.
 
virtual const ElementGraphics< F > * graphics () const
 
uint & tag ()
 Returns the tag.
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 

Detailed Description

template<class F = Real>
class hp2Dedge::Edge< F >

An edge of a 2D FEM edge element quad.

Useful for integration over the edge, in a bilinearform or linearform. E.g. hp2Dedge::Dirichlet need it to calculate the local coefficients.

Definition at line 162 of file quad.hh.

Member Typedef Documentation

◆ type

template<class F >
typedef F concepts::Element< F >::type
inherited

Definition at line 54 of file element.hh.

Constructor & Destructor Documentation

◆ Edge()

template<class F = Real>
hp2Dedge::Edge< F >::Edge ( const Quad< F > &  elm,
const ushort  k 
)

Constructor

Parameters
elmelement, on which the edge lies
knumber of edge, 0 - left, 1 - upper, 2 - right, 3 - bottom

Member Function Documentation

◆ chi()

template<class F = Real>
concepts::Real2d hp2Dedge::Edge< F >::chi ( const Real  t) const
inline

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

Definition at line 232 of file quad.hh.

◆ diffElement()

template<class F = Real>
Real hp2Dedge::Edge< F >::diffElement ( const Real  t) const
inline

Computes the differential element for integration over [-1,1].

Definition at line 260 of file quad.hh.

◆ direction()

template<class F = Real>
const ushort hp2Dedge::Edge< F >::direction ( ) const
inline

Returns direction of edge on reference quad [0,1]^2, 0 - x, 1 - y.

Definition at line 197 of file quad.hh.

◆ edge()

template<class F = Real>
const ushort hp2Dedge::Edge< F >::edge ( ) const
inline

Returns number of the edge.

Definition at line 185 of file quad.hh.

◆ elm()

template<class F = Real>
const Quad< F > & hp2Dedge::Edge< F >::elm ( ) const
inline

Returns element.

Definition at line 191 of file quad.hh.

◆ graphics()

◆ info()

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

Returns information in an output stream.

Reimplemented from concepts::OutputOperator.

◆ integration()

template<class F = Real>
const concepts::QuadratureRule1d * hp2Dedge::Edge< F >::integration ( ) const
inline

Returns the integration rule.

Definition at line 216 of file quad.hh.

◆ jacobian()

template<class F = Real>
concepts::MapReal2d hp2Dedge::Edge< F >::jacobian ( const Real  t) const
inline

Computes the Jacobian matrix of element transformation on the edge

Definition at line 241 of file quad.hh.

◆ jacobianDeterminant()

template<class F = Real>
Real hp2Dedge::Edge< F >::jacobianDeterminant ( const Real  t) const
inline

Computes the determinant of the Jacobian.

Definition at line 254 of file quad.hh.

◆ jacobianInverse()

template<class F = Real>
concepts::MapReal2d hp2Dedge::Edge< F >::jacobianInverse ( const Real  t) const
inline

Computes the inverse of the Jacobian.

Definition at line 248 of file quad.hh.

◆ localCoords()

template<class F = Real>
concepts::Real2d hp2Dedge::Edge< F >::localCoords ( const Real  t) const
inline

coordinate of point on the edge inside reference element [0,1]^2

Definition at line 222 of file quad.hh.

◆ shpfct()

template<class F = Real>
const hp2D::KarniadakisDeriv2 * hp2Dedge::Edge< F >::shpfct ( ) const
inline

Definition at line 211 of file quad.hh.

◆ sign()

template<class F = Real>
const Real hp2Dedge::Edge< F >::sign ( ) const
inline

Returns sign of outer normal vector, e.g. left edge -1, right edge +1.

Definition at line 203 of file quad.hh.

◆ support()

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

Definition at line 179 of file quad.hh.

◆ T()

template<class F = Real>
virtual const concepts::TMatrixBase< F > & hp2Dedge::Edge< F >::T ( ) const
inlinevirtual

T-Matrix of the appropiate Quad, not used.

Implements concepts::Element< F >.

Definition at line 268 of file quad.hh.

◆ tag()

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

Returns the tag.

Definition at line 66 of file element.hh.

◆ vertex()

template<class F = Real>
concepts::Real2d hp2Dedge::Edge< F >::vertex ( uint  i) const
inline

Definition at line 173 of file quad.hh.


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