Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)

Public Types | Public Member Functions | Static Public Member Functions | Protected Attributes | Static Protected Attributes
hp2D::IntegrableQuad Class Reference

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

#include <quad.hh>

Inheritance diagram for hp2D::IntegrableQuad:
Inheritance graph
[legend]
Collaboration diagram for hp2D::IntegrableQuad:
Collaboration graph
[legend]

List of all members.

Public Types

enum  intFormType { ZERO, ONE, TWO, THREE }
 Integration form, which determines terms coming from integration over reference element. More...

Public Member Functions

concepts::Real2d chi (const Real x, const Real y) const
 Computes the element map, the map from the reference element to the real geometry of the quad.
const concepts::Quad2dSubdivisiongetStrategy () const
 Returns the subdivision strategy of the underlying cell of this element.
 IntegrableQuad (concepts::Quad2d &cell)
const concepts::QuadratureRuleintegrationX () const
 Returns the integration rule in x direction.
const concepts::QuadratureRuleintegrationY () const
 Returns the integration rule in y direction.
concepts::MapReal2d jacobian (const Real x, const Real y) const
 Computes the Jacobian.
Real jacobianDeterminant (const Real x, const Real y) const
 Computes the determinant of the Jacobian.
concepts::MapReal2d jacobianInverse (const Real x, const Real y) const
 Computes the inverse of the Jacobian.
virtual bool quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const
 Delivers a quadrature point.
void setStrategy (const concepts::Quad2dSubdivision *strategy=0) throw (concepts::StrategyChange)
 Sets the subdivision strategy of the underlying cell of this element.

Static Public Member Functions

static concepts::QuadRuleFactoryrule ()
 Access to the quadrature rule, which is valid for all elements of this type (hp2D::IntegrableQuad).

Protected Attributes

concepts::Quad2dcell_
 The cell.
std::auto_ptr
< concepts::QuadratureRule
intX_
 The integration rules.
std::auto_ptr
< concepts::QuadratureRule
intY_

Static Protected Attributes

static concepts::QuadRuleFactory rule_

Detailed Description

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

Definition at line 29 of file quad.hh.


Member Enumeration Documentation

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

Enumerator:
ZERO 
ONE 
TWO 
THREE 

Definition at line 27 of file integral.hh.


Constructor & Destructor Documentation

hp2D::IntegrableQuad::IntegrableQuad ( concepts::Quad2d cell)

Member Function Documentation

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 55 of file quad.hh.

Here is the call graph for this function:

const concepts::Quad2dSubdivision* hp2D::IntegrableQuad::getStrategy ( ) const [inline]

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

Definition at line 96 of file quad.hh.

Here is the call graph for this function:

const concepts::QuadratureRule* hp2D::IntegrableQuad::integrationX ( ) const [inline]

Returns the integration rule in x direction.

Definition at line 34 of file quad.hh.

const concepts::QuadratureRule* hp2D::IntegrableQuad::integrationY ( ) const [inline]

Returns the integration rule in y direction.

Definition at line 38 of file quad.hh.

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

Computes the Jacobian.

Definition at line 60 of file quad.hh.

Here is the call graph for this function:

Real hp2D::IntegrableQuad::jacobianDeterminant ( const Real  x,
const Real  y 
) const [inline]

Computes the determinant of the Jacobian.

Definition at line 75 of file quad.hh.

Here is the call graph for this function:

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

Computes the inverse of the Jacobian.

Definition at line 68 of file quad.hh.

Here is the call graph for this function:

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.

static concepts::QuadRuleFactory& hp2D::IntegrableQuad::rule ( ) [inline, static]

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, exactDtN.cc, inhomDirichletBCs.cc, inhomDirichletBCsLagrange.cc, inhomNeumannBCs.cc, and RobinBCs.cc.

Definition at line 48 of file quad.hh.

void hp2D::IntegrableQuad::setStrategy ( const concepts::Quad2dSubdivision strategy = 0) throw (concepts::StrategyChange) [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 87 of file quad.hh.

Here is the call graph for this function:


Member Data Documentation

The cell.

Definition at line 110 of file quad.hh.

The integration rules.

Definition at line 106 of file quad.hh.

Definition at line 107 of file quad.hh.

Definition at line 113 of file quad.hh.


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

Home | Doxygen Documentation | Tutorials | Developer Tools (restricted)