DUNE-FEM (unstable)

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

constructor More...

#include <dune/fem/quadrature/cachingpointlist.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 QuadraturePointIterator< ThisIteratorType
 type of iterator
 
typedef ElementIntegrationPointList< GridPartType, codimension, IntegrationTraits > NonConformingQuadratureType
 type of quadrature used for non-conforming intersections
 
enum  Side
 inside and outside flags
 
typedef IntegrationTraits::IntegrationPointListType IntegrationPointListType
 type of the integration point list
 

Public Member Functions

 CachingPointList (const GridPartType &gridPart, const IntersectionType &intersection, const QuadratureKeyType &quadKey, const typename Base ::Side side)
 constructor More...
 
const CoordinateTypepoint (const size_t i) const
 obtain coordinates of i-th integration point More...
 
int twistId () const
 returns the twistId, i.e. [0,...,7] More...
 
size_t cachingPoint (const size_t quadraturePoint) const
 map quadrature points to caching points More...
 
size_t interpolationPoint (const size_t quadraturePoint) const
 map quadrature points to interpolation points More...
 
bool isInterpolationQuadrature (const size_t numShapeFunctions) const
 check if quadrature is interpolation quadrature 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...
 

Static Public Member Functions

static constexpr bool twisted ()
 returns true if cachingPoint is not the identity mapping More...
 

Protected Member Functions

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

Detailed Description

template<typename GridPartImp, class IntegrationTraits>
class Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >

constructor

Note
The CachingPointList requires the grid part to get twist information for TwistUtility (see also ElementIntegrationPointList<GridPartImp,1>).
Parameters
[in]gridPartgrid partition
[in]intersectionintersection
[in]quadKeydesired order of the quadrature or other means of quadrature identification
[in]sideeither INSIDE or OUTSIDE; codim-0 entity for which the ElementQuadrature shall be created

Constructor & Destructor Documentation

◆ CachingPointList()

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

constructor

Note
The CachingPointList requires the grid part to get twist information for TwistUtility (see also ElementIntegrationPointList<GridPartImp,1>).
Parameters
[in]gridPartgrid partition
[in]intersectionintersection
[in]quadKeydesired order of the quadrature or other means of quadrature identification
[in]sideeither INSIDE or OUTSIDE; codim-0 entity for which the ElementQuadrature shall be created

Member Function Documentation

◆ cachingPoint()

template<typename GridPartImp , class IntegrationTraits >
size_t Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::cachingPoint ( const size_t  quadraturePoint) const
inline

map quadrature points to caching points

For codim-1 entites, the mapping consists of two stages:

  • Consider the twist to get the quadrature point number on the face of the (codim-0) reference element,
  • Map the twisted quadrature point number to the caching point number.
Parameters
[in]quadraturePointnumber of quadrature point to map to a caching point

◆ 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

◆ interpolationPoint()

template<typename GridPartImp , class IntegrationTraits >
size_t Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::interpolationPoint ( const size_t  quadraturePoint) const
inline

map quadrature points to interpolation points

Parameters
[in]quadraturePointnumber of quadrature point to map to an interpolation point

◆ isInterpolationQuadrature()

template<typename GridPartImp , class IntegrationTraits >
bool Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::isInterpolationQuadrature ( const size_t  numShapeFunctions) const
inline

check if quadrature is interpolation quadrature

Parameters
[in]numShapeFunctionsnumber of shapeFunctions that has to match number of quadrature points or number of internal interpolation points

◆ 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<typename GridPartImp , class IntegrationTraits >
const CoordinateType& Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::point ( const 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

◆ twisted()

template<typename GridPartImp , class IntegrationTraits >
static constexpr bool Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::twisted ( )
inlinestaticconstexpr

returns true if cachingPoint is not the identity mapping

◆ twistId()

template<typename GridPartImp , class IntegrationTraits >
int Dune::Fem::CachingPointList< GridPartImp, 1, IntegrationTraits >::twistId ( ) const
inline

returns the twistId, i.e. [0,...,7]


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)