DUNE-FEM (unstable)

Dune::Fem::QuadratureRulesFactory< FieldImp, dim > Class Template Reference

quadrature implementation based on the standard DUNE quadratures More...

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

Public Types

enum  
 to be revised, look at caching quad
 

Public Member Functions

 QuadratureRulesFactory (const GeometryType &geometry, const int order, const size_t id)
 constructor filling the list of points and weights More...
 
int order () const
 obtain order of the integration point list More...
 
GeometryType geometryType () const
 
const FieldType & weight (size_t i) const
 obtain weight of i-th integration point More...
 
const CoordinateType & point (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 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 Member Functions

static unsigned int maxOrder ()
 maximal order of available quadratures
 

Protected Member Functions

void addQuadraturePoint (const CoordinateType &point, const FieldType weight)
 Adds a point-weight pair to the quadrature. More...
 
void setIntegrationPoints (std::vector< CoordinateType > &&points)
 Overwrites integration point list

 

Detailed Description

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

quadrature implementation based on the standard DUNE quadratures

Though a factory by name, this is a quadrature implementation using the standard quadratures from DUNE grid to generate a list of quadrature points.

Constructor & Destructor Documentation

◆ QuadratureRulesFactory()

template<typename FieldImp , int dim>
Dune::Fem::QuadratureRulesFactory< FieldImp, dim >::QuadratureRulesFactory ( const GeometryType geometry,
const int  order,
const size_t  id 
)
inline

constructor filling the list of points and weights

Parameters
[in]geometrygeometry type for which a quadrature is desired
[in]orderdesired order (provided by the user)
[in]idunique identifier (provided by QuadratureProvider)

References Dune::Fem::QuadratureImp< FieldImp, dim >::addQuadraturePoint(), Dune::QuadratureType::GaussLobatto, Dune::Fem::QuadratureRulesFactory< FieldImp, dim >::order(), and Dune::QuadratureRule< ct, dim >::order().

Member Function Documentation

◆ addQuadraturePoint()

template<typename FieldImp , int dim>
void Dune::Fem::QuadratureImp< FieldImp, dim >::addQuadraturePoint ( const CoordinateType point,
const FieldType  weight 
)
inlineprotectedinherited

Adds a point-weight pair to the quadrature.

This method allows derived classes to add quadrature points (and their respective weights) to the list. This mehtod should only be used within the constructor of the derived class.

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

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

◆ geometryType()

template<typename FieldImp , int dim>
GeometryType Dune::Fem::QuadratureRulesFactory< FieldImp, dim >::geometryType ( ) const
inlinevirtual

◆ id()

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

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

◆ nop()

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

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>
int Dune::Fem::QuadratureRulesFactory< FieldImp, dim >::order ( ) const
inlinevirtual

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

Implements Dune::Fem::IntegrationPointListImp< FieldImp, dim >.

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

◆ point()

template<typename FieldImp , int dim>
const CoordinateType& Dune::Fem::IntegrationPointListImp< FieldImp, dim >::point ( size_t  i) const
inlineinherited

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().

◆ weight()

template<typename FieldImp , int dim>
const FieldType& Dune::Fem::QuadratureImp< FieldImp, dim >::weight ( size_t  i) const
inlineinherited

obtain weight of i-th integration point

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

Referenced by 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)