1#ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH
2#define DUNE_FEM_ELEMENTPOINTLISTBASE_HH
4#include <dune/geometry/referenceelements.hh>
17 template<
class Gr
idPartImp,
int codim,
class IntegrationTraits >
18 class ElementPointListBase;
21 template<
class Gr
idPartImp,
class IntegrationTraits >
40 static const int dimension = GridPartType::dimension;
107 return quadImp().geometryType();
114 return quadImp().geometryType();
121 return quadImp().geometryType();
145 return quadraturePoint;
151 return quadraturePoint;
155 inline bool twisted ()
const {
return false; }
187 template<
class Gr
idPartImp,
int codim,
class IntegrationTraits >
224 const GeometryType &faceGeo,
227 : quad_( faceGeo, quadKey ),
228 elementGeometry_( elementGeo ),
242 elementGeometry_( elementGeo ),
265 return quad_.point( i );
304 return elementGeometry_;
309 return quadraturePoint;
314 return quadraturePoint;
318 inline bool twisted ()
const {
return false; }
328 return localFaceIndex_;
349 static const bool isCube =
353 static const bool isSimplex =
357 if( isCube || isSimplex )
362 return Dune::GeometryTypes::cube(
dimension-1 );
367 return Dune::GeometryTypes::simplex(
dimension-1 );
370 else if( elementGeo.isNone() )
374 if( elementGeo.dim() == 2 )
376 return Dune::GeometryTypes::cube( 1 );
380 return Dune::GeometryTypes::none( elementGeo.dim()-1 );
385 assert( ! elementGeo.isNone() );
386 typedef Dune::ReferenceElements< RealType, dimension > RefElements;
387 return RefElements::general( elementGeo ).type( face,
codimension );
393 GeometryType elementGeometry_;
Definition: bindguard.hh:11
specialize with 'true' for if the codimension 0 entity of the grid part has only one possible geometr...
Definition: gridpart/common/capabilities.hh:29
ElementPointListBase.
Definition: elementpointlistbase.hh:189
IntegrationPointListType::CoordinateType LocalCoordinateType
Definition: elementpointlistbase.hh:212
ElementPointListBase(const GeometryType &elementGeo, const GeometryType &faceGeo, const int localFaceIndex, const QuadratureKeyType &quadKey)
constructor
Definition: elementpointlistbase.hh:223
Side
inside and outside flags
Definition: elementpointlistbase.hh:197
@ OUTSIDE
Definition: elementpointlistbase.hh:197
@ INSIDE
Definition: elementpointlistbase.hh:197
IntegrationTraits::CoordinateType CoordinateType
Definition: elementpointlistbase.hh:211
bool twisted() const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:318
size_t localCachingPoint(const size_t quadraturePoint) const
Definition: elementpointlistbase.hh:312
GridPartImp GridPartType
type of the grid partition
Definition: elementpointlistbase.hh:194
int twistId() const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:321
const LocalCoordinateType & localPoint(size_t i) const
obtain local coordinates of i-th integration point
Definition: elementpointlistbase.hh:263
IntegrationPointListType::QuadratureKeyType QuadratureKeyType
Definition: elementpointlistbase.hh:214
const IntegrationPointListType & quadImp() const
obtain the actual implementation of the quadrature
Definition: elementpointlistbase.hh:338
ElementPointListBase(const GeometryType &elementGeo, const int localFaceIndex, const QuadratureKeyType &quadKey)
constructor
Definition: elementpointlistbase.hh:238
static const int dimension
dimension of the grid
Definition: elementpointlistbase.hh:206
size_t cachingPoint(const size_t quadraturePoint) const
Definition: elementpointlistbase.hh:307
static GeometryType getFaceGeometry(const GeometryType &elementGeo, const int face)
Definition: elementpointlistbase.hh:344
size_t nop() const
obtain the number of integration points
Definition: elementpointlistbase.hh:248
int localFaceIndex() const
Definition: elementpointlistbase.hh:326
IntegrationTraits::IntegrationPointListType IntegrationPointListType
type of the integration point list
Definition: elementpointlistbase.hh:209
GeometryType geometry() const
obtain GeometryType for this integration point list
Definition: elementpointlistbase.hh:284
int cachingPointStart() const
Definition: elementpointlistbase.hh:324
int nCachingPoints() const
Definition: elementpointlistbase.hh:323
static const int codimension
codimension of the element integration point list
Definition: elementpointlistbase.hh:200
size_t id() const
obtain the identifier of the integration point list
Definition: elementpointlistbase.hh:270
GridPartType::ctype RealType
coordinate type
Definition: elementpointlistbase.hh:203
int order() const
obtain order of the integration point list
Definition: elementpointlistbase.hh:277
GeometryType elementGeometry() const
obtain GeometryType of the corresponding codim-0 the integration point list belongs to
Definition: elementpointlistbase.hh:302
Definition: elementpointlistbase.hh:23
Side
inside and outside flags
Definition: elementpointlistbase.hh:31
@ INSIDE
Definition: elementpointlistbase.hh:31
ElementPointListBase(const IntegrationPointListType &ipList)
constructor
Definition: elementpointlistbase.hh:64
GridPartType::ctype RealType
coordinate type
Definition: elementpointlistbase.hh:37
int nCachingPoints() const
Definition: elementpointlistbase.hh:165
size_t cachingPoint(const size_t quadraturePoint) const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:143
int twistId() const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:158
GeometryType elementGeometry() const
obtain GeometryType of the corresponding codim-0 the integration point list belongs to
Definition: elementpointlistbase.hh:137
GeometryType type() const
Definition: elementpointlistbase.hh:112
bool twisted() const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:155
GeometryType geometry() const
Definition: elementpointlistbase.hh:105
IntegrationPointListType quad_
Definition: elementpointlistbase.hh:181
int order() const
obtain order of the integration point list
Definition: elementpointlistbase.hh:98
const LocalCoordinateType & localPoint(size_t i) const
obtain local coordinates of i-th integration point
Definition: elementpointlistbase.hh:84
ElementPointListBase(const GeometryType &geometry, const QuadratureKeyType &quadKey)
constructor
Definition: elementpointlistbase.hh:55
IntegrationTraits::IntegrationPointListType IntegrationPointListType
type of the integration point list
Definition: elementpointlistbase.hh:43
IntegrationPointListType::CoordinateType LocalCoordinateType
Definition: elementpointlistbase.hh:46
int localFaceIndex() const
Definition: elementpointlistbase.hh:160
size_t id() const
obtain the identifier of the integration point list
Definition: elementpointlistbase.hh:91
const IntegrationPointListType & quadImp() const
obtain the actual implementation of the quadrature
Definition: elementpointlistbase.hh:175
size_t localCachingPoint(const size_t quadraturePoint) const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:149
GeometryType geometryType() const
Definition: elementpointlistbase.hh:119
GridPartImp GridPartType
type of the grid partition
Definition: elementpointlistbase.hh:28
IntegrationPointListType::QuadratureKeyType QuadratureKeyType
Definition: elementpointlistbase.hh:48
int cachingPointStart() const
Definition: elementpointlistbase.hh:166
IntegrationTraits::CoordinateType CoordinateType
Definition: elementpointlistbase.hh:45
size_t nop() const
obtain the number of integration points
Definition: elementpointlistbase.hh:69