Class documentation of Concepts

Loading...
Searching...
No Matches

#include <quad.hh>

Inheritance diagram for hp2D::IntegrableQuad:
concepts::IntegrationCell hp2D::BaseQuad< Real > hp2D::BaseQuad< H > hp2D::BaseQuad< F > hp2D::Quad< Real > hp2D::Quad< H > hp2D::Quad< F > hp2Dedge::Quad< F > hp3D::NeumannTraceElement3d< F >

Public Types

enum  intFormType { ZERO , ONE , TWO , THREE }
 

Public Member Functions

 IntegrableQuad (concepts::QuadNd &cell)
 
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
 

Static Public Member Functions

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

Protected Attributes

std::unique_ptr< concepts::QuadratureRule2dintRule_
 The integration rules.
 
concepts::QuadNdcell_
 The cell.
 

Static Protected Attributes

static std::unique_ptr< concepts::QuadRuleFactoryTensor2dfactory_
 

Detailed Description

Class holding the quadrature rule and the cell of a quadrilateral element.

Definition at line 31 of file quad.hh.

Member Enumeration Documentation

◆ intFormType

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

Definition at line 29 of file integral.hh.

Member Function Documentation

◆ chi()

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

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.

◆ factory()

static std::unique_ptr< concepts::QuadRuleFactoryTensor2d > & hp2D::IntegrableQuad::factory ( )
inlinestatic

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 ( )
inlinestatic

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
inline

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
inline

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.

◆ intRule()

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

Definition at line 135 of file quad.hh.

◆ inverseLaplace()

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

Definition at line 66 of file quad.hh.

◆ jacobian()

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

Computes the Jacobian.

Definition at line 45 of file quad.hh.

◆ jacobianInverse()

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

Computes the inverse of the Jacobian.

Definition at line 59 of file quad.hh.

◆ quadraturePoint()

virtual bool hp2D::IntegrableQuad::quadraturePoint ( uint  i,
intPoint p,
intFormType  form = ZERO,
bool  localCoord = false 
) const
virtual

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.

◆ setStrategy()

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

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.

Member Data Documentation

◆ cell_

concepts::QuadNd& hp2D::IntegrableQuad::cell_
protected

The cell.

Definition at line 161 of file quad.hh.

◆ factory_

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

Definition at line 165 of file quad.hh.

◆ intRule_

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

The integration rules.

Definition at line 158 of file quad.hh.


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