Class documentation of Concepts

Loading...
Searching...
No Matches

#include <quad.hh>

Inheritance diagram for hp2D::Quad< F >:
hp2D::BaseQuad< F > hp2D::QuadShapeFunctions hp2D::Element< F > hp2D::IntegrableQuad concepts::ElementWithCell< F > concepts::IntegrationCell concepts::Element< F > concepts::OutputOperator hp3D::NeumannTraceElement3d< F >

Public Types

typedef F type
 
enum  intFormType { ZERO , ONE , TWO , THREE }
 

Public Member Functions

 Quad (concepts::QuadNd &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
 
virtual const concepts::ElementGraphics< F > * graphics () const
 Returns element graphics class.
 
void recomputeShapefunctions ()
 
void recomputeShapefunctions (const uint nq[2])
 
void edgeP (const uint i, const uint &p)
 Set polynomial degree of edge i to p.
 
const uint edgeP (const uint i) const
 
virtual const concepts::Quadsupport () const
 
virtual concepts::Real3d vertex (uint i) const
 
virtual const concepts::QuadNdcell () 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 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::Quad2dSubdivisiongetStrategy () const
 
virtual bool quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const
 
const concepts::QuadratureRule2d *const intRule () const
 
const ushort * p () const
 
const concepts::Karniadakis< 1, 0 > * shpfctX () const
 Returns the shape functions in x direction.
 
const concepts::Karniadakis< 1, 0 > * shpfctY () const
 Returns the shape functions in y direction.
 
const concepts::Karniadakis< 1, 1 > * shpfctDX () const
 Returns the derivatives of the shape functions in x direction.
 
const concepts::Karniadakis< 1, 1 > * shpfctDY () const
 Returns the shape functions in y direction.
 

Static Public Member Functions

static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & factory ()
 
static concepts::QuadRuleFactoryTensor2dfactory_rp ()
 

Protected Member Functions

virtual std::ostream & info (std::ostream &os) const
 Returns information in an output stream.
 
void computeShapefunctions_ (const concepts::QuadratureRule2d *intRule)
 gets the shapefunctions, used in both constructors
 

Protected Attributes

concepts::TMatrix< F > T_
 T matrix of the element.
 
std::unique_ptr< concepts::QuadratureRule2dintRule_
 The integration rules.
 
concepts::QuadNdcell_
 The cell.
 

Static Protected Attributes

static std::unique_ptr< concepts::QuadRuleFactoryTensor2dfactory_
 

Detailed Description

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

A 2D FEM element: a quad.

The reference shape functions are products of the polynomials of Karniadakis and Sherwin. The index of the shape functions rises first over the polynomials in local x-direction.

Definition at line 271 of file quad.hh.

Member Typedef Documentation

◆ type

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

Definition at line 81 of file element.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

◆ Quad()

template<class F = Real>
hp2D::Quad< F >::Quad ( concepts::QuadNd cell,
const ushort *  p,
concepts::TColumn< F > *  T0,
concepts::TColumn< F > *  T1 
)

Constructor

Parameters
cellCell on which the element is defined
pPolynomial degree (might be anisotropic)
T0Part of the T matrix
T1Part of the T matrix

Member Function Documentation

◆ appendT()

template<class F >
void hp2D::Element< F >::appendT ( concepts::TColumn< F > *  T)
inlineinherited

Appends the T columns to the T matrix.

Definition at line 43 of file element.hh.

◆ cell()

template<class F = Real>
virtual const concepts::QuadNd & hp2D::BaseQuad< F >::cell ( ) const
inlinevirtualinherited

Returns the cell on which the element is built.

Implements hp2D::Element< F >.

Definition at line 196 of file quad.hh.

◆ chi()

concepts::Real2d hp2D::IntegrableQuad::chi ( const Real  x,
const Real  y 
) const
inlineinherited

Computes the element map, the map from the reference element to the real geometry of the quad.

The reference element is the unit quad.

Definition at line 40 of file quad.hh.

◆ computeShapefunctions_()

void hp2D::QuadShapeFunctions::computeShapefunctions_ ( const concepts::QuadratureRule2d intRule)
protectedinherited

gets the shapefunctions, used in both constructors

gets the shapefunctions, used in both constructors

◆ edgeP() [1/2]

template<class F = Real>
const uint hp2D::Quad< F >::edgeP ( const uint  i) const
inline

Definition at line 295 of file quad.hh.

◆ edgeP() [2/2]

template<class F = Real>
void hp2D::Quad< F >::edgeP ( const uint  i,
const uint &  p 
)
inline

Set polynomial degree of edge i to p.

Definition at line 292 of file quad.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.

◆ factory()

static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & hp2D::IntegrableQuad::factory ( )
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.

Examples
BGT_0.cc, RobinBCs.cc, elasticity2D_tutorial.cc, exactDtN.cc, howToGetStarted.cc, inhomDirichletBCs.cc, inhomDirichletBCsLagrange.cc, inhomNeumannBCs.cc, and parallelizationTutorial.cc.

Definition at line 144 of file quad.hh.

◆ factory_rp()

static concepts::QuadRuleFactoryTensor2d * hp2D::IntegrableQuad::factory_rp ( )
inlinestaticinherited

Access to the quadrature rule ( as factory() ) through a raw pointer.

Definition at line 148 of file quad.hh.

◆ getStrategy()

const concepts::Quad2dSubdivision * hp2D::IntegrableQuad::getStrategy ( ) const
inlineinherited

Returns the subdivision strategy of the underlying cell of this element.

Definition at line 121 of file quad.hh.

◆ gramDeterminantRoot()

Real hp2D::IntegrableQuad::gramDeterminantRoot ( const Real  x,
const Real  y 
) const
inlineinherited

Computes the sqaure root of the Gram determinant multiplied by 0.25 to incorporate the transformation factor for the integration weights.

Definition at line 91 of file quad.hh.

◆ graphics()

template<class F = Real>
virtual const concepts::ElementGraphics< F > * hp2D::Quad< F >::graphics ( ) const
virtual

Returns element graphics class.

Implements hp2D::BaseQuad< F >.

◆ info()

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

Returns information in an output stream.

Reimplemented from hp2D::BaseQuad< F >.

Reimplemented in hp3D::NeumannTraceElement3d< F >.

◆ intRule()

const concepts::QuadratureRule2d *const hp2D::IntegrableQuad::intRule ( ) const
inlineinherited

Definition at line 135 of file quad.hh.

◆ inverseLaplace()

concepts::MapReal2d hp2D::IntegrableQuad::inverseLaplace ( const Real  x,
const Real  y 
) const
inlineinherited

Definition at line 66 of file quad.hh.

◆ jacobian()

concepts::MapReal2d hp2D::IntegrableQuad::jacobian ( const Real  x,
const Real  y 
) const
inlineinherited

Computes the Jacobian.

Definition at line 45 of file quad.hh.

◆ jacobianInverse()

concepts::MapReal2d hp2D::IntegrableQuad::jacobianInverse ( const Real  x,
const Real  y 
) const
inlineinherited

Computes the inverse of the Jacobian.

Definition at line 59 of file quad.hh.

◆ p()

const ushort * hp2D::QuadShapeFunctions::p ( ) const
inlineinherited

Returns the polynomial degree. The returned array has 2 elements.

Definition at line 229 of file quad.hh.

◆ quadraturePoint()

virtual bool hp2D::IntegrableQuad::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::Quad< F >::recomputeShapefunctions ( )

Recompute shape functions, e.g. for other abscissas redefined through setIntegrationRule

◆ setStrategy()

void hp2D::IntegrableQuad::setStrategy ( const concepts::Quad2dSubdivision strategy = 0)
inlineinherited

Sets the subdivision strategy of the underlying cell of this element. It calls Quad2d::setStrategy.

Parameters
strategyPointer to an instance of a subdivision strategy.
Exceptions
StrategyChangeif the change is not allowed (the change is not allowed if there are children present)

Definition at line 109 of file quad.hh.

◆ shpfctDX()

const concepts::Karniadakis< 1, 1 > * hp2D::QuadShapeFunctions::shpfctDX ( ) const
inlineinherited

Returns the derivatives of the shape functions in x direction.

Definition at line 240 of file quad.hh.

◆ shpfctDY()

const concepts::Karniadakis< 1, 1 > * hp2D::QuadShapeFunctions::shpfctDY ( ) const
inlineinherited

Returns the shape functions in y direction.

Definition at line 243 of file quad.hh.

◆ shpfctX()

const concepts::Karniadakis< 1, 0 > * hp2D::QuadShapeFunctions::shpfctX ( ) const
inlineinherited

Returns the shape functions in x direction.

Definition at line 232 of file quad.hh.

◆ shpfctY()

const concepts::Karniadakis< 1, 0 > * hp2D::QuadShapeFunctions::shpfctY ( ) const
inlineinherited

Returns the shape functions in y direction.

Definition at line 236 of file quad.hh.

◆ support()

template<class F = Real>
virtual const concepts::Quad & hp2D::BaseQuad< F >::support ( ) const
inlinevirtualinherited

Returns the topological support of the element. Possible supports for an element are quadrilaterals and triangles.

Implements hp2D::Element< F >.

Definition at line 194 of file quad.hh.

◆ T()

template<class F >
virtual const concepts::TMatrix< F > & hp2D::Element< F >::T ( ) const
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.

◆ 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>
virtual concepts::Real3d hp2D::BaseQuad< F >::vertex ( uint  i) const
inlinevirtualinherited

Returns the coordinates of the ith vertex of this element.

Parameters
iIndex of the vertex

Implements hp2D::Element< F >.

Definition at line 195 of file quad.hh.

Member Data Documentation

◆ cell_

concepts::QuadNd& hp2D::IntegrableQuad::cell_
protectedinherited

The cell.

Definition at line 161 of file quad.hh.

◆ factory_

std::unique_ptr<concepts::QuadRuleFactoryTensor2d> hp2D::IntegrableQuad::factory_
staticprotectedinherited

Definition at line 165 of file quad.hh.

◆ intRule_

std::unique_ptr<concepts::QuadratureRule2d> hp2D::IntegrableQuad::intRule_
protectedinherited

The integration rules.

Definition at line 158 of file quad.hh.

◆ T_

template<class F >
concepts::TMatrix<F> hp2D::Element< F >::T_
protectedinherited

T matrix of the element.

Definition at line 57 of file element.hh.


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