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

Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Attributes | Static Private Attributes
hp2D::Quad< F > Class Template Reference

A 2D FEM element: a quad. More...

#include <quad.hh>

Inheritance diagram for hp2D::Quad< F >:
Inheritance graph
[legend]
Collaboration diagram for hp2D::Quad< F >:
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...
typedef F type

Public Member Functions

void appendT (concepts::TColumn< F > *T)
 Appends the T columns to the T matrix.
virtual const concepts::Quad2dcell () const
 Returns the cell on which the element is built.
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.
void edgeP (const uint i, const uint &p)
 Set polynomial degree of edge i to p.
const uint edgeP (const uint i) const
Real3d elemMap (const Real coord_local) const
Real3d elemMap (const Real2d &coord_local) const
Real3d elemMap (const Real3d &coord_local) const
const concepts::Quad2dSubdivisiongetStrategy () const
 Returns the subdivision strategy of the underlying cell of this element.
virtual const
concepts::ElementGraphics< F > * 
graphics () const
 Returns element graphics class.
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.
const ushort * p () const
 Returns the polynomial degree.
 Quad (concepts::Quad2d &cell, const ushort *p, concepts::TColumn< F > *T0, concepts::TColumn< F > *T1)
 Constructor.
virtual bool quadraturePoint (uint i, intPoint &p, intFormType form=ZERO, bool localCoord=false) const
 Delivers a quadrature point.
void recomputeShapefunctions ()
 Recompute shape functions, e.g.
void recomputeShapefunctions (const uint nq[2])
void setStrategy (const concepts::Quad2dSubdivision *strategy=0) throw (concepts::StrategyChange)
 Sets the subdivision strategy of the underlying cell of this element.
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.
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.
virtual const concepts::Quadsupport () const
 Returns the topological support of the element.
virtual const
concepts::TMatrix< F > & 
T () const
 Returns the T matrix of the element.
uint & tag ()
 Returns the tag.
virtual concepts::Real3d vertex (uint i) const
 Returns the coordinates of the ith vertex of this element.
virtual ~Quad ()

Static Public Member Functions

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

Protected Member Functions

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

Protected Attributes

concepts::Quad2dcell_
 The cell.
std::auto_ptr
< concepts::QuadratureRule
intX_
 The integration rules.
std::auto_ptr
< concepts::QuadratureRule
intY_
concepts::TMatrix< F > T_
 T matrix of the element.

Static Protected Attributes

static concepts::QuadRuleFactory rule_

Private Attributes

const uint * edges_ [4]
 Polynomial degree of edges.

Static Private Attributes

static std::auto_ptr
< concepts::ElementGraphics< F > > 
graphics_
 Appropiate element graphics object.

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


Member Typedef Documentation

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

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

template<class F = Real>
hp2D::Quad< F >::Quad ( concepts::Quad2d 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
template<class F = Real>
virtual hp2D::Quad< F >::~Quad ( ) [virtual]

Member Function Documentation

template<class F>
void hp2D::Element< F >::appendT ( concepts::TColumn< F > *  T) [inline, inherited]

Appends the T columns to the T matrix.

Definition at line 41 of file element.hh.

template<class F = Real>
virtual const concepts::Quad2d& hp2D::BaseQuad< F >::cell ( ) const [inline, virtual, inherited]

Returns the cell on which the element is built.

Possible are tetrahedrons, hexahedron, prims and pyramids.

Implements hp2D::Element< F >.

Definition at line 135 of file quad.hh.

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

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:

void hp2D::QuadShapeFunctions::computeShapefunctions_ ( const concepts::QuadratureRule intX,
const concepts::QuadratureRule intY 
) [protected, inherited]

gets the shapefunctions, used in both constructors

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

Definition at line 234 of file quad.hh.

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

template<typename F>
Real3d concepts::ElementWithCell< F >::elemMap ( const Real2d coord_local) const [inline, inherited]

Definition at line 87 of file element.hh.

template<typename F>
Real3d concepts::ElementWithCell< F >::elemMap ( const Real  coord_local) const [inline, inherited]

Definition at line 83 of file element.hh.

template<typename F>
Real3d concepts::ElementWithCell< F >::elemMap ( const Real3d coord_local) const [inline, inherited]

Definition at line 91 of file element.hh.

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

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:

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

Returns element graphics class.

Implements hp2D::BaseQuad< F >.

template<class F = Real>
virtual std::ostream& hp2D::Quad< F >::info ( std::ostream &  os) const [protected, virtual]

Returns information in an output stream.

Reimplemented from hp2D::BaseQuad< F >.

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

Returns the integration rule in x direction.

Definition at line 34 of file quad.hh.

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

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, inherited]

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, inherited]

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, inherited]

Computes the inverse of the Jacobian.

Definition at line 68 of file quad.hh.

Here is the call graph for this function:

const ushort* hp2D::QuadShapeFunctions::p ( ) const [inline, inherited]

Returns the polynomial degree.

The returned array has 2 elements.

Definition at line 172 of file quad.hh.

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

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.

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

Recompute shape functions, e.g.

for other abscissas redefined through setIntegrationRule

template<class F = Real>
void hp2D::Quad< F >::recomputeShapefunctions ( const uint  nq[2])
static concepts::QuadRuleFactory& hp2D::IntegrableQuad::rule ( ) [inline, static, inherited]

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, inherited]

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:

const concepts::Karniadakis<1,1>* hp2D::QuadShapeFunctions::shpfctDX ( ) const [inline, inherited]

Returns the derivatives of the shape functions in x direction.

Definition at line 183 of file quad.hh.

const concepts::Karniadakis<1,1>* hp2D::QuadShapeFunctions::shpfctDY ( ) const [inline, inherited]

Returns the shape functions in y direction.

Definition at line 186 of file quad.hh.

const concepts::Karniadakis<1,0>* hp2D::QuadShapeFunctions::shpfctX ( ) const [inline, inherited]

Returns the shape functions in x direction.

Definition at line 175 of file quad.hh.

const concepts::Karniadakis<1,0>* hp2D::QuadShapeFunctions::shpfctY ( ) const [inline, inherited]

Returns the shape functions in y direction.

Definition at line 179 of file quad.hh.

template<class F = Real>
virtual const concepts::Quad& hp2D::BaseQuad< F >::support ( ) const [inline, virtual, inherited]

Returns the topological support of the element.

Possible supports for an element are quadrilaterals and triangles.

Implements hp2D::Element< F >.

Definition at line 133 of file quad.hh.

template<class F>
virtual const concepts::TMatrix<F>& hp2D::Element< F >::T ( ) const [inline, virtual, inherited]

Returns the T matrix of the element.

Implements concepts::Element< F >.

Definition at line 38 of file element.hh.

template<class F>
uint& concepts::Element< F >::tag ( ) [inline, inherited]

Returns the tag.

Definition at line 65 of file element.hh.

template<class F = Real>
virtual concepts::Real3d hp2D::BaseQuad< F >::vertex ( uint  i) const [inline, virtual, inherited]

Returns the coordinates of the ith vertex of this element.

Parameters:
iIndex of the vertex

Implements hp2D::Element< F >.

Definition at line 134 of file quad.hh.


Member Data Documentation

The cell.

Definition at line 110 of file quad.hh.

template<class F = Real>
const uint* hp2D::Quad< F >::edges_[4] [private]

Polynomial degree of edges.

Definition at line 245 of file quad.hh.

template<class F = Real>
std::auto_ptr<concepts::ElementGraphics<F> > hp2D::Quad< F >::graphics_ [static, private]

Appropiate element graphics object.

Definition at line 243 of file quad.hh.

std::auto_ptr<concepts::QuadratureRule> hp2D::IntegrableQuad::intX_ [protected, inherited]

The integration rules.

Definition at line 106 of file quad.hh.

std::auto_ptr<concepts::QuadratureRule> hp2D::IntegrableQuad::intY_ [protected, inherited]

Definition at line 107 of file quad.hh.

Definition at line 113 of file quad.hh.

template<class F>
concepts::TMatrix<F> hp2D::Element< F >::T_ [protected, inherited]

T matrix of the element.

Definition at line 57 of file element.hh.


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

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