1 #ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH 2 #define DUNE_FEM_ELEMENTPOINTLISTBASE_HH 4 #include <dune/geometry/referenceelements.hh> 15 template<
class Gr
idPartImp,
int codim,
class IntegrationTraits >
19 template<
class Gr
idPartImp,
class IntegrationTraits >
38 static const int dimension = GridPartType::dimension;
52 : quad_( geometry, order )
94 return quadImp().geometryType();
101 return quadImp().geometryType();
108 return quadImp().geometryType();
131 return quadraturePoint;
141 const IntegrationPointListType &
quadImp ()
const 147 IntegrationPointListType quad_;
153 template<
class Gr
idPartImp,
int codim,
class IntegrationTraits >
191 const GeometryType &faceGeo,
194 : quad_( faceGeo, order ),
195 elementGeometry_( elementGeo ),
196 localFaceIndex_( localFaceIndex )
209 elementGeometry_( elementGeo ),
210 localFaceIndex_( localFaceIndex )
213 const QuadraturePointWrapperType
operator[] (
size_t i )
const 237 return quad_.point( i );
276 return elementGeometry_;
281 return quadraturePoint;
291 const IntegrationPointListType &
quadImp()
const 298 return localFaceIndex_;
304 typedef Dune::ReferenceElements< RealType, dimension > RefElements;
305 return RefElements::general( elementGeo ).type( face, codimension );
309 IntegrationPointListType quad_;
310 GeometryType elementGeometry_;
318 #endif // #ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH IntegrationPointListType::CoordinateType LocalCoordinateType
Definition: elementpointlistbase.hh:44
Definition: elementpointlistbase.hh:20
int order() const
obtain order of the integration point list
Definition: elementpointlistbase.hh:249
int localFaceIndex() const
Definition: elementpointlistbase.hh:296
size_t cachingPoint(const size_t quadraturePoint) const
Definition: elementpointlistbase.hh:279
static GeometryType getFaceGeometry(const GeometryType &elementGeo, const int face)
Definition: elementpointlistbase.hh:302
Side
inside and outside flags
Definition: elementpointlistbase.hh:163
int order() const
obtain order of the integration point list
Definition: elementpointlistbase.hh:85
const IntegrationPointListType & quadImp() const
obtain the actual implementation of the quadrature
Definition: elementpointlistbase.hh:291
size_t nop() const
obtain the number of integration points
Definition: elementpointlistbase.hh:220
ElementPointListBase(const GeometryType &elementGeo, const int localFaceIndex, const int order)
constructor
Definition: elementpointlistbase.hh:205
ElementPointListBase(const GeometryType &elementGeo, const GeometryType &faceGeo, const int localFaceIndex, const int order)
constructor
Definition: elementpointlistbase.hh:190
wrapper for a (Quadrature,int) pair
Definition: quadrature.hh:40
const QuadraturePointWrapperType operator[](size_t i) const
Definition: elementpointlistbase.hh:213
ElementPointListBase.
Definition: elementpointlistbase.hh:16
GeometryType geometry() const
obtain GeometryType for this integration point list
Definition: elementpointlistbase.hh:256
size_t id() const
obtain the identifier of the integration point list
Definition: elementpointlistbase.hh:78
const IntegrationPointListType & quadImp() const
obtain the actual implementation of the quadrature
Definition: elementpointlistbase.hh:141
IntegrationPointListType::CoordinateType LocalCoordinateType
Definition: elementpointlistbase.hh:178
GridPartImp GridPartType
type of the grid partition
Definition: elementpointlistbase.hh:26
size_t id() const
obtain the identifier of the integration point list
Definition: elementpointlistbase.hh:242
size_t nop() const
obtain the number of integration points
Definition: elementpointlistbase.hh:56
GridPartImp GridPartType
type of the grid partition
Definition: elementpointlistbase.hh:160
IntegrationTraits::CoordinateType CoordinateType
Definition: elementpointlistbase.hh:177
Definition: coordinate.hh:4
GeometryType elementGeometry() const
obtain GeometryType of the corresponding codim-0 the integration point list belongs to ...
Definition: elementpointlistbase.hh:274
GeometryType elementGeometry() const
obtain GeometryType of the corresponding codim-0 the integration point list belongs to ...
Definition: elementpointlistbase.hh:124
GridPartType::ctype RealType
coordinate type
Definition: elementpointlistbase.hh:35
Definition: elementpointlistbase.hh:163
GridPartType::ctype RealType
coordinate type
Definition: elementpointlistbase.hh:169
ElementPointListBase(const GeometryType &geometry, int order)
constructor
Definition: elementpointlistbase.hh:51
const LocalCoordinateType & localPoint(size_t i) const
obtain local coordinates of i-th integration point
Definition: elementpointlistbase.hh:235
IntegrationTraits::IntegrationPointListType IntegrationPointListType
type of the integration point list
Definition: elementpointlistbase.hh:175
GeometryType geometryType() const
Definition: elementpointlistbase.hh:106
size_t cachingPoint(const size_t quadraturePoint) const
Definition: elementpointlistbase.hh:129
static const int dimension
dimension of the grid
Definition: elementpointlistbase.hh:172
IntegrationTraits::IntegrationPointListType IntegrationPointListType
type of the integration point list
Definition: elementpointlistbase.hh:41
GeometryType type() const
Definition: elementpointlistbase.hh:99
GeometryType geometry() const
Definition: elementpointlistbase.hh:92
static const int codimension
codimension of the element integration point list
Definition: elementpointlistbase.hh:166
QuadraturePointWrapper< This > QuadraturePointWrapperType
the type of the quadrature point
Definition: elementpointlistbase.hh:181
Side
inside and outside flags
Definition: elementpointlistbase.hh:29
const LocalCoordinateType & localPoint(size_t i) const
obtain local coordinates of i-th integration point
Definition: elementpointlistbase.hh:71
IntegrationTraits::CoordinateType CoordinateType
Definition: elementpointlistbase.hh:43
Definition: elementpointlistbase.hh:163