DUNE-FEM (unstable)
#include <dune/fem/quadrature/lumpingquadrature.hh>
Public Types | |
enum | |
to be revised, look at caching quad | |
Public Member Functions | |
LumpingQuadrature (const GeometryType >, int order, int id) | |
constructor filling the list of points and weights. More... | |
virtual GeometryType | geometryType () const |
virtual int | order () const |
obtain order of the integration point list More... | |
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< ElementCoordinateType > | interpolationPoints (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 std::size_t | 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
class Dune::Fem::LumpingQuadrature< FieldImp, geometryId >
Define a lumping quadrature for all geometries. Note, however, that this may not make sense for anything else than simplices or maybe hexagonal grids. For simplicial meshes the quadrature formula is exact on linear polynomials and hence the quadrature error is quadratic in the mesh-size. A mass-matrix assembled with the caching quadrature will be diagonal in the context of Lagrange spaces. Generally, it is a bad idea to use mass-lumping for anything else than linear (or maybe bilinear) finite elements.
There are probably much more efficient ways to perform mass-lumping. This "quadrature" approach is convenient, however, as it can simply be plugged into existing code by replacing the quadrature rules.
Constructor & Destructor Documentation
◆ LumpingQuadrature()
|
inline |
constructor filling the list of points and weights.
- Parameters
-
[in] gt geometry type for which a quadrature is desired [in] order order (between 1 and 6) [in] id unique identifier, ignored
References Dune::Fem::QuadratureImp< FieldImp, Dune::GeometryType(geometryId).dim()>::addQuadraturePoint(), Dune::Geo::ReferenceElements< ctype_, dim >::general(), Dune::FloatCmp::gt(), Dune::GeometryType::id(), and Dune::Fem::QuadratureImp< FieldImp, Dune::GeometryType(geometryId).dim()>::weight().
Member Function Documentation
◆ addQuadraturePoint()
|
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.
◆ geometryType()
|
inlinevirtual |
◆ id()
|
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()
|
inlineinherited |
obtain the number of integration points
- Returns
- number of integration points within this list
Referenced by Dune::Fem::IntegrationPointListImp< FieldImp, dim >::point().
◆ order()
|
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 >.
◆ point()
|
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] i number 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()
|
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] i number of the integration point, 0 <= i < nop()
- Returns
- weight of the i-th integration point
The documentation for this class was generated from the following file:
- dune/fem/quadrature/lumpingquadrature.hh