DUNE-FEM (unstable)

Dune::Fem::IntegrationPointListImp< FieldImp, dim > Class Template Referenceabstract

Generic implementation of an IntegrationPointList. More...

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

Public Types

enum  
 to be revised, look at caching quad
 
typedef FieldImp FieldType
 field type
 
typedef FieldVector< FieldType, dim > CoordinateType
 type of local coordinates
 

Public Member Functions

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...
 
size_t id () const
 obtain the identifier of the integration point list More...
 
virtual int order () const =0
 obtain order of the integration point list More...
 
virtual GeometryType geometryType () const =0
 obtain GeometryType for this integration point list More...
 
virtual std::vector< ElementCoordinateTypeinterpolationPoints (const int reqDim) const
 returns list of element interpolation points for a given face quadrature
 
virtual bool isFaceInterpolationQuadrature (const size_t numShapeFunctions) const
 return true if quadrature is also a set of interpolation points for a given number of shape functions
 

Static Public Attributes

static const int dimension = dim
 dimension of quadrature
 

Protected Member Functions

 IntegrationPointListImp (size_t id)
 Constructor. More...
 
void addIntegrationPoint (const CoordinateType &point)
 Adds an integration point to the list. More...
 
void setIntegrationPoints (std::vector< CoordinateType > &&points)
 Overwrites integration point list

 

Detailed Description

template<typename FieldImp, int dim>
class Dune::Fem::IntegrationPointListImp< FieldImp, dim >

Generic implementation of an IntegrationPointList.

An integration point list is simply a list of points, given in local coordinates, i.e., coordinates within the reference element.

Note
Integration point lists do not change over time. It can safely be assumed that they always return the same points in the same order.

Constructor & Destructor Documentation

◆ IntegrationPointListImp()

template<typename FieldImp , int dim>
Dune::Fem::IntegrationPointListImp< FieldImp, dim >::IntegrationPointListImp ( size_t  id)
inlineexplicitprotected

Constructor.

The constructor simply creates an empty point list and stores the specified identifier.

Note
The constructors of derived classes should fill the integration point list via addIntegrationPoint.
The identifier of an integration point list must be globally unique. Even integration point lists for different dimensions must have different identifiers.
Parameters
[in]idunique identifier of the integration point list (provided by QuadratureProvider)

Member Function Documentation

◆ addIntegrationPoint()

template<typename FieldImp , int dim>
void Dune::Fem::IntegrationPointListImp< FieldImp, dim >::addIntegrationPoint ( const CoordinateType point)
inlineprotected

Adds an integration point to the list.

This method allows derived classes to add integration points to the list. This mehtod should only be used within the constructor of the derived class.

References Dune::Fem::IntegrationPointListImp< FieldImp, dim >::point().

◆ geometryType()

template<typename FieldImp , int dim>
virtual GeometryType Dune::Fem::IntegrationPointListImp< FieldImp, dim >::geometryType ( ) const
pure virtual

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.

Returns
GeometryType for this integration point list

Implemented in Dune::Fem::LagrangePointListInterface< FieldImp, dim, maxPolOrder >, Dune::Fem::LumpingQuadrature< FieldImp, geometryId >, and Dune::Fem::QuadratureRulesFactory< FieldImp, dim >.

Referenced by Dune::Fem::TwistMapperCreator< ct, dim >::TwistMapperCreator().

◆ id()

template<typename FieldImp , int dim>
size_t Dune::Fem::IntegrationPointListImp< FieldImp, dim >::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

Referenced by Dune::Fem::TwistProvider< ct, dim >::getTwistStorage().

◆ nop()

template<typename FieldImp , int dim>
size_t Dune::Fem::IntegrationPointListImp< FieldImp, dim >::nop ( ) const
inline

obtain the number of integration points

Returns
number of integration points within this list

Referenced by Dune::Fem::IntegrationPointListImp< FieldImp, dim >::point().

◆ order()

template<typename FieldImp , int dim>
virtual int Dune::Fem::IntegrationPointListImp< FieldImp, dim >::order ( ) const
pure virtual

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.

Returns
the order of the integration point list

Implemented in Dune::Fem::LagrangePointListInterface< FieldImp, dim, maxPolOrder >, Dune::Fem::LumpingQuadrature< FieldImp, geometryId >, and Dune::Fem::QuadratureRulesFactory< FieldImp, dim >.

◆ point()

template<typename FieldImp , int dim>
const CoordinateType& Dune::Fem::IntegrationPointListImp< FieldImp, dim >::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::IntegrationPointListImp< FieldImp, dim >::nop().

Referenced by Dune::Fem::IntegrationPointListImp< FieldImp, dim >::addIntegrationPoint(), and Dune::Fem::QuadratureImp< FieldImp, dim >::addQuadraturePoint().


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 4, 22:30, 2024)