DUNE-FEM (unstable)

Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature > Class Template Reference

actual interface class for integration point lists More...

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

Public Types

enum  
 to be revised, look at caching quad
 
typedef Traits::IntegrationPointListType IntegrationPointListType
 type of integration point list implementation
 
typedef IntegrationPointListType::CoordinateType CoordinateType
 type of coordinate
 

Public Member Functions

template<class FactoryTraits >
 IntegrationPointList (const FactoryTraits traits, const GeometryType &geometryType, const typename FactoryTraits::QuadratureKeyType &quadKey)
 create a quadrature for a given geometry type and order More...
 
template<class FactoryTraits >
 IntegrationPointList (const FactoryTraits traits, const GeometryType &geometryType, const GeometryType &elementGeometry, const typename FactoryTraits::QuadratureKeyType &quadKey)
 create a quadrature for a given geometry type and order More...
 
 IntegrationPointList (const IntegrationPointListType &ipList)
 create an integration point list from an implementation More...
 
 IntegrationPointList (const IntegrationPointListStorageType &ipListPtr)
 create an integration point list from an implementation More...
 
 IntegrationPointList (const IntegrationPointList &org)=default
 copy constructor More...
 
const IntegrationPointListTypeipList () const
 obtain a reference the actual implementation More...
 
int nop () const
 obtain the number of integration points More...
 
const CoordinateTypepoint (size_t i) const
 obtain 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 geometryType () const
 obtain GeometryType for this integration point list More...
 
auto interpolationPoints (const int reqDim) const
 returns list of element interpolation points for a given face quadrature
 
bool isFaceInterpolationQuadrature (const size_t numShapeFunctions) const
 return true if quadrature is also a set of interpolation points for the given shape functions
 
WeightReturnType weight (size_t i) const
 obtain weight of i-th integration point (if quadrature, else 1.0) More...
 

Detailed Description

template<typename FieldImp, int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
class Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >

actual interface class for integration point lists

IntegrationPointList is a proxy for the actual implementations of the integration point lists. During construction, the IntegrationPointList object is configured with an appropriate implementation object from the QuadratureProvider (monostate pattern).

The design goal is minimization of construction time. The actual implementation can be created once and reused whenever it is needed. Moreover, this layout insulates the user from all initialization and storage stuff.

Note
The difference between integration point lists and quadratures is that quadratures have weights.

Constructor & Destructor Documentation

◆ IntegrationPointList() [1/5]

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
template<class FactoryTraits >
Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::IntegrationPointList ( const FactoryTraits  traits,
const GeometryType geometryType,
const typename FactoryTraits::QuadratureKeyType &  quadKey 
)
inline

create a quadrature for a given geometry type and order

This constructor creates a quadrature for the specified geometry which is capable of integrating polynomials up the given order exactly.

Note
The order of the quadrature may be higher than the requested one.
Parameters
[in]geometryTypegeometry type of the requested quadrature
[in]orderorder of the requested quadrature

◆ IntegrationPointList() [2/5]

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
template<class FactoryTraits >
Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::IntegrationPointList ( const FactoryTraits  traits,
const GeometryType geometryType,
const GeometryType elementGeometry,
const typename FactoryTraits::QuadratureKeyType &  quadKey 
)
inline

create a quadrature for a given geometry type and order

This constructor creates a quadrature for the specified geometry which is capable of integrating polynomials up the given order exactly.

Note
The order of the quadrature may be higher than the requested one.
Parameters
[in]geometryTypegeometry type of the requested quadrature
[in]elementGeometrygeometry type of element that resulting quadrature is used for (in case of face quadratures)
[in]orderorder of the requested quadrature

◆ IntegrationPointList() [3/5]

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::IntegrationPointList ( const IntegrationPointListType ipList)
inline

create an integration point list from an implementation

This constructor creates an integration point list from a given implementation.

Note
This constructor is provided mainly for testing purposes.
Parameters
[in]ipListimplementation of the integration point list

◆ IntegrationPointList() [4/5]

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::IntegrationPointList ( const IntegrationPointListStorageType &  ipListPtr)
inline

create an integration point list from an implementation

This constructor creates an integration point list from a given implementation.

Note
This constructor is provided mainly for agglomeration quadratures
Parameters
[in]ipListPtrimplementation of the integration point list

◆ IntegrationPointList() [5/5]

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::IntegrationPointList ( const IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature > &  org)
inlinedefault

copy constructor

Parameters
[in]orgintegration point list to be copied

Member Function Documentation

◆ geometryType()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
GeometryType Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::geometryType ( ) const
inline

obtain GeometryType for this integration point list

Integration point lists are specified in local coordinates, i.e., coordinates with respect to the reference element. Hence, each integration point list is only valid for one type of geometry, i.e., for one reference element. The type can be retrieved via this method.

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

References Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::ipList().

◆ id()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
size_t Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::id ( ) const
inline

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

References Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::ipList().

◆ ipList()

◆ nop()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
int Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::nop ( ) const
inline

◆ order()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
int Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::order ( ) const
inline

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

References Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::ipList().

◆ point()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
const CoordinateType & Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::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

References Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::ipList().

◆ weight()

template<typename FieldImp , int dim, template< class, int > class IntegrationTraits, bool isQuadrature = false>
WeightReturnType Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::weight ( size_t  i) const
inline

obtain weight of i-th integration point (if quadrature, else 1.0)

This method returns the weight of the i-th integration point for 0 <= i < nop() within the quadrature.

Note
The integration point can be obtained via the point() method.
The quadrature weights sum up to the volume of the reference element.
Parameters
[in]inumber of the integration point, 0 <= i < nop()
Returns
weight of the i-th integration point

References Dune::Fem::IntegrationPointList< FieldImp, dim, IntegrationTraits, isQuadrature >::ipList().


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.111.3 (Nov 20, 23:30, 2024)