DUNE-FEM (unstable)

Dune::Fem::ElementIntegrationPointList< GridPartImp, 1, IntegrationTraits > Class Template Reference

constructor More...

#include <dune/fem/quadrature/elementpointlist.hh>

Public Types

typedef GridPartImp GridPartType
 type of the grid partition
 
typedef Base::CoordinateType CoordinateType
 Type of coordinates in codim-0 reference element.
 
typedef Base::QuadratureKeyType QuadratureKeyType
 type of quadrature identifier on user side (default is the order of quadrature)
 
typedef GridPartType::IntersectionIteratorType IntersectionIteratorType
 Type of the intersection iterator.
 
typedef QuadraturePointWrapper< ThisQuadraturePointWrapperType
 type of the quadrature point
 
typedef QuadraturePointIterator< ThisIteratorType
 type of iterator
 
typedef This NonConformingQuadratureType
 type quadrature for use on non-conforming intersections
 
enum  Side
 inside and outside flags
 
typedef GridPartType::ctype RealType
 coordinate type
 
typedef IntegrationTraits::IntegrationPointListType IntegrationPointListType
 type of the integration point list
 

Public Member Functions

 ElementIntegrationPointList (const GridPartType &gridPart, const IntersectionType &intersection, const QuadratureKeyType &quadKey, const typename Base ::Side side)
 constructor More...
 
const CoordinateTypepoint (size_t i) const
 obtain coordinates of i-th integration point More...
 
size_t nop () const
 obtain the number of integration points More...
 
const LocalCoordinateType & localPoint (size_t i) const
 obtain local coordinates of i-th integration point More...
 
size_t id () const
 obtain the identifier of the integration point list More...
 
int order () const
 obtain order of the integration point list More...
 
GeometryType geometry () const
 obtain GeometryType for this integration point list
 
GeometryType elementGeometry () const
 obtain GeometryType of the corresponding codim-0 the integration point list belongs to More...
 
int twistId () const
 convenience implementation for Dune::Fem::CachingInterface
 

Static Public Member Functions

static constexpr bool twisted ()
 convenience implementation for Dune::Fem::CachingInterface
 

Static Public Attributes

static const int codimension
 codimension of the element integration point list
 

Protected Member Functions

const IntegrationPointListTypequadImp () const
 obtain the actual implementation of the quadrature More...
 

Detailed Description

template<class GridPartImp, class IntegrationTraits>
class Dune::Fem::ElementIntegrationPointList< GridPartImp, 1, IntegrationTraits >

constructor

Parameters
[in]gridPartgrid partition (a dummy here)
[in]intersectionintersection
[in]quadKeyquadrature key, i.e. desired order of the quadrature
[in]sideeither INSIDE or OUTSIDE; codim-0 entity for which the ElementQuadrature shall be created
Note
This code assumes that the codim-0 entity is either a simplex or a cube (otherwise elementGeometry() returns a wrong geometry).

Constructor & Destructor Documentation

◆ ElementIntegrationPointList()

template<class GridPartImp , class IntegrationTraits >
Dune::Fem::ElementIntegrationPointList< GridPartImp, 1, IntegrationTraits >::ElementIntegrationPointList ( const GridPartType gridPart,
const IntersectionType &  intersection,
const QuadratureKeyType quadKey,
const typename Base ::Side  side 
)
inline

constructor

Parameters
[in]gridPartgrid partition (a dummy here)
[in]intersectionintersection
[in]quadKeyquadrature key, i.e. desired order of the quadrature
[in]sideeither INSIDE or OUTSIDE; codim-0 entity for which the ElementQuadrature shall be created
Note
This code assumes that the codim-0 entity is either a simplex or a cube (otherwise elementGeometry() returns a wrong geometry).

Member Function Documentation

◆ elementGeometry()

GeometryType Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::elementGeometry ( ) const
inlineinherited

obtain GeometryType of the corresponding codim-0 the integration point list belongs to

An element integration point list can return the coordinates of integration points with resepct to the codim-0 reference element and the reference element corresponding to the subentity the quadrature actually lives on. This method returns the geometry of the codim-0 entity.

Note
Calling this method yields a virtual function call, so do not call this method unnecessarily.
Returns
GeometryType for this integration point list

◆ id()

size_t Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::id ( ) const
inlineinherited

obtain the identifier of the integration point list

The identifier of an integration point list must be globally unique. Even integration point lists for different dimensions must have different identifiers.

Note
Quadratures are considered distinct if they differ in one of the following points: geometry type, order, dimension or implementation.
Returns
globally unique identifier of the integration point list

◆ localPoint()

const LocalCoordinateType& Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::localPoint ( size_t  i) const
inlineinherited

obtain local coordinates of i-th integration point

This method returns a reference to the local coordinates of the i-th integration point for 0 <= i < nop(). Here, local coordinates means coordinates with respect to the reference element of the subentity.

Parameters
[in]inumber of the integration point, 0 <= i < nop()
Returns
reference to i-th integration point

◆ nop()

size_t Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::nop ( ) const
inlineinherited

obtain the number of integration points

Returns
number of integration points within this list

◆ order()

int Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::order ( ) const
inlineinherited

obtain order of the integration point list

The order of a quadrature is the maximal polynomial degree that is guaranteed to be integrated exactly by the quadrature.

In case of an integration point list, the definition of this value is left to the implementor.

Note
Calling this method yields a virtual function call, so do not call this method unnecessarily.
Returns
the order of the integration point list

◆ point()

template<class GridPartImp , class IntegrationTraits >
const CoordinateType& Dune::Fem::ElementIntegrationPointList< GridPartImp, 1, IntegrationTraits >::point ( size_t  i) const
inline

obtain coordinates of i-th integration point

This method returns a reference to the coordinates of the i-th integration point for 0 <= i < nop(). The integration point is given in local coordinates, i.e., coordinates with respect to the reference element.

Parameters
[in]inumber of the integration point, 0 <= i < nop()
Returns
reference to i-th integration point

◆ quadImp()

const IntegrationPointListType& Dune::Fem::ElementPointListBase< GridPartImp, codim, IntegrationTraits >::quadImp ( ) const
inlineprotectedinherited

obtain the actual implementation of the quadrature

Note
This method may only be used in derived classes.
Returns
a reference to the actual implementation of the quadrature

The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 16, 22:29, 2024)