dune-fem  2.4.1-rc
elementpointlistbase.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH
2 #define DUNE_FEM_ELEMENTPOINTLISTBASE_HH
3 
4 #include <dune/geometry/referenceelements.hh>
5 
7 
8 namespace Dune
9 {
10 
11  namespace Fem
12  {
13 
15  template< class GridPartImp, int codim, class IntegrationTraits >
17 
18 
19  template< class GridPartImp, class IntegrationTraits >
20  class ElementPointListBase< GridPartImp, 0, IntegrationTraits >
21  {
23 
24  public:
26  typedef GridPartImp GridPartType;
27 
29  enum Side { INSIDE, OUTSIDE };
30 
32  static const int codimension = 0;
33 
35  typedef typename GridPartType::ctype RealType;
36 
38  static const int dimension = GridPartType::dimension;
39 
41  typedef typename IntegrationTraits::IntegrationPointListType IntegrationPointListType;
42 
43  typedef typename IntegrationTraits::CoordinateType CoordinateType;
44  typedef typename IntegrationPointListType::CoordinateType LocalCoordinateType;
45 
51  ElementPointListBase ( const GeometryType &geometry, int order )
52  : quad_( geometry, order )
53  {}
54 
56  size_t nop () const
57  {
58  return quadImp().nop();
59  }
60 
71  const LocalCoordinateType &localPoint( size_t i ) const
72  {
73  return quadImp().point( i );
74  }
75 
78  size_t id () const
79  {
80  return quadImp().id();
81  }
82 
85  int order () const
86  {
87  return quadImp().order();
88  }
89 
92  GeometryType geometry () const
93  {
94  return quadImp().geometryType();
95  }
96 
99  GeometryType type () const
100  {
101  return quadImp().geometryType();
102  }
103 
106  GeometryType geometryType () const
107  {
108  return quadImp().geometryType();
109  }
110 
124  GeometryType elementGeometry () const
125  {
126  return quadImp().geometry();
127  }
128 
129  size_t cachingPoint( const size_t quadraturePoint ) const
130  {
131  return quadraturePoint;
132  }
133 
134  protected:
141  const IntegrationPointListType &quadImp () const
142  {
143  return quad_;
144  }
145 
146  private:
147  IntegrationPointListType quad_;
148  };
149 
150 
151 
153  template< class GridPartImp, int codim, class IntegrationTraits >
155  {
157 
158  public:
160  typedef GridPartImp GridPartType;
161 
163  enum Side { INSIDE, OUTSIDE };
164 
166  static const int codimension = codim;
167 
169  typedef typename GridPartType::ctype RealType;
170 
172  static const int dimension = GridPartType::dimension;
173 
175  typedef typename IntegrationTraits::IntegrationPointListType IntegrationPointListType;
176 
177  typedef typename IntegrationTraits::CoordinateType CoordinateType;
178  typedef typename IntegrationPointListType::CoordinateType LocalCoordinateType;
179 
182 
190  ElementPointListBase ( const GeometryType &elementGeo,
191  const GeometryType &faceGeo,
192  const int localFaceIndex,
193  const int order )
194  : quad_( faceGeo, order ),
195  elementGeometry_( elementGeo ),
196  localFaceIndex_( localFaceIndex )
197  {}
198 
205  ElementPointListBase ( const GeometryType &elementGeo,
206  const int localFaceIndex,
207  const int order )
208  : quad_( getFaceGeometry( elementGeo, localFaceIndex ), order ),
209  elementGeometry_( elementGeo ),
210  localFaceIndex_( localFaceIndex )
211  {}
212 
213  const QuadraturePointWrapperType operator[] ( size_t i ) const
214  {
215  return QuadraturePointWrapperType( *this, i );
216  }
217 
220  size_t nop () const
221  {
222  return quadImp().nop();
223  }
224 
235  const LocalCoordinateType &localPoint ( size_t i ) const
236  {
237  return quad_.point( i );
238  }
239 
242  size_t id () const
243  {
244  return quadImp().id();
245  }
246 
249  int order () const
250  {
251  return quadImp().order();
252  }
253 
256  GeometryType geometry () const
257  {
258  return quadImp().geo();
259  }
260 
274  GeometryType elementGeometry () const
275  {
276  return elementGeometry_;
277  }
278 
279  size_t cachingPoint( const size_t quadraturePoint ) const
280  {
281  return quadraturePoint;
282  }
283 
284  protected:
291  const IntegrationPointListType &quadImp() const
292  {
293  return quad_;
294  }
295 
296  int localFaceIndex () const
297  {
298  return localFaceIndex_;
299  }
300 
301  static GeometryType
302  getFaceGeometry ( const GeometryType &elementGeo, const int face )
303  {
304  typedef Dune::ReferenceElements< RealType, dimension > RefElements;
305  return RefElements::general( elementGeo ).type( face, codimension );
306  }
307 
308  private:
309  IntegrationPointListType quad_;
310  GeometryType elementGeometry_;
311  int localFaceIndex_;
312  };
313 
314  } // namespace Fem
315 
316 } // namespace Dune
317 
318 #endif // #ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH
IntegrationPointListType::CoordinateType LocalCoordinateType
Definition: elementpointlistbase.hh:44
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