DUNE-FEM (unstable)

elementpointlistbase.hh
1#ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH
2#define DUNE_FEM_ELEMENTPOINTLISTBASE_HH
3
4#include <dune/geometry/referenceelements.hh>
5
6#include <dune/fem/quadrature/quadrature.hh>
7
8#include <dune/fem/gridpart/common/capabilities.hh>
9
10namespace Dune
11{
12
13 namespace Fem
14 {
15
17 template< class GridPartImp, int codim, class IntegrationTraits >
18 class ElementPointListBase;
19
20
21 template< class GridPartImp, class IntegrationTraits >
22 class ElementPointListBase< GridPartImp, 0, IntegrationTraits >
23 {
25
26 public:
28 typedef GridPartImp GridPartType;
29
31 enum Side { INSIDE, OUTSIDE };
32
34 static const int codimension = 0;
35
37 typedef typename GridPartType::ctype RealType;
38
40 static const int dimension = GridPartType::dimension;
41
43 typedef typename GridPartType::template Codim< 0 >::EntityType EntityType;
44
46 typedef typename IntegrationTraits::IntegrationPointListType IntegrationPointListType;
47
48 typedef typename IntegrationTraits::CoordinateType CoordinateType;
49 typedef typename IntegrationPointListType::CoordinateType LocalCoordinateType;
50
51 //typedef typename IntegrationPointListType :: QuadratureKeyType QuadratureKeyType;
52
58 template <class QuadratureKey>
59 ElementPointListBase ( const GeometryType &geometry, const QuadratureKey& quadKey )
60 : quad_( geometry, quadKey )
61 {}
62
69 : quad_( ipList )
70 {}
71
73 size_t nop () const
74 {
75 return quadImp().nop();
76 }
77
88 const LocalCoordinateType &localPoint( size_t i ) const
89 {
90 return quadImp().point( i );
91 }
92
95 size_t id () const
96 {
97 return quadImp().id();
98 }
99
102 int order () const
103 {
104 return quadImp().order();
105 }
106
109 GeometryType geometry () const
110 {
111 return quadImp().geometryType();
112 }
113
116 GeometryType type () const
117 {
118 return quadImp().geometryType();
119 }
120
124 {
125 return quadImp().geometryType();
126 }
127
142 {
143 return quadImp().geometry();
144 }
145
147 size_t cachingPoint( const size_t quadraturePoint ) const
148 {
149 return quadraturePoint;
150 }
151
153 size_t localCachingPoint( const size_t quadraturePoint ) const
154 {
155 return quadraturePoint;
156 }
157
159 static constexpr bool twisted () { return false; }
160
162 inline int twistId () const { return 0; }
163
164 int localFaceIndex () const
165 {
166 return 0;
167 }
168
169 inline int nCachingPoints () const { return nop(); }
170 inline int cachingPointStart () const { return 0; }
171
172 protected:
179 const IntegrationPointListType &quadImp () const
180 {
181 return quad_;
182 }
183
184 protected:
186 };
187
188
189
191 template< class GridPartImp, int codim, class IntegrationTraits >
193 {
195
196 public:
198 typedef GridPartImp GridPartType;
199
201 enum Side { INSIDE, OUTSIDE };
202
204 static const int codimension = codim;
205
207 typedef typename GridPartType::ctype RealType;
208
210 static const int dimension = GridPartType::dimension;
211
213 typedef typename IntegrationTraits::IntegrationPointListType IntegrationPointListType;
214
215 typedef typename IntegrationTraits::CoordinateType CoordinateType;
216 typedef typename IntegrationPointListType::CoordinateType LocalCoordinateType;
217
225 template <class QuadratureKeyType>
227 const GeometryType &faceGeo,
228 const int localFaceIndex,
229 const QuadratureKeyType& quadKey )
230 : quad_( faceGeo, quadKey ),
231 elementGeometry_( elementGeo ),
232 localFaceIndex_( localFaceIndex )
233 {}
234
241 template <class QuadratureKeyType>
243 const int localFaceIndex,
244 const QuadratureKeyType& quadKey )
245 : quad_( getFaceGeometry( elementGeo, localFaceIndex ), quadKey ),
246 elementGeometry_( elementGeo ),
247 localFaceIndex_( localFaceIndex )
248 {}
249
252 size_t nop () const
253 {
254 return quadImp().nop();
255 }
256
267 const LocalCoordinateType &localPoint ( size_t i ) const
268 {
269 return quad_.point( i );
270 }
271
274 size_t id () const
275 {
276 return quadImp().id();
277 }
278
281 int order () const
282 {
283 return quadImp().order();
284 }
285
289 {
290 return quadImp().geo();
291 }
292
307 {
308 return elementGeometry_;
309 }
310
311 size_t cachingPoint( const size_t quadraturePoint ) const
312 {
313 return quadraturePoint;
314 }
315
316 size_t localCachingPoint( const size_t quadraturePoint ) const
317 {
318 return quadraturePoint;
319 }
320
322 static constexpr bool twisted () { return false; }
323
325 inline int twistId () const { return 0; }
326
327 inline int nCachingPoints () const { return nop(); }
328 inline int cachingPointStart () const { return 0; }
329
330 int localFaceIndex () const
331 {
332 return localFaceIndex_;
333 }
334
335 protected:
343 {
344 return quad_;
345 }
346
347 static GeometryType
348 getFaceGeometry ( const GeometryType &elementGeo, const int face )
349 {
350 // for cube and simplex geom types the dim-1 geom type
351 // is also cube or simplex
352
353 static const bool isCube =
356
357 static const bool isSimplex =
360
361 if( isCube || isSimplex )
362 {
363 assert( elementGeo.dim() == dimension );
364 if( isCube )
365 {
367 }
368 else
369 {
370 assert( isSimplex );
372 }
373 }
374 else if( elementGeo.isNone() )
375 {
376 // if cell geometry is none and dim is 2 then the
377 // face is a normal edge which is of type cube
378 if( elementGeo.dim() == 2 )
379 {
380 return Dune::GeometryTypes::cube( 1 );
381 }
382 else
383 {
384 return Dune::GeometryTypes::none( elementGeo.dim()-1 );
385 }
386 }
387 else // use reference element to determine type
388 {
389 assert( ! elementGeo.isNone() );
391 return RefElements::general( elementGeo ).type( face, codimension );
392 }
393 }
394
395 private:
397 GeometryType elementGeometry_;
398 int localFaceIndex_;
399 };
400
401 } // namespace Fem
402
403} // namespace Dune
404
405#endif // #ifndef DUNE_FEM_ELEMENTPOINTLISTBASE_HH
Side
inside and outside flags
Definition: elementpointlistbase.hh:201
static constexpr bool twisted()
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:322
GridPartImp GridPartType
type of the grid partition
Definition: elementpointlistbase.hh:198
int twistId() const
convenience implementation for Dune::Fem::CachingInterface
Definition: elementpointlistbase.hh:325
ElementPointListBase(const GeometryType &elementGeo, const int localFaceIndex, const QuadratureKeyType &quadKey)
constructor
Definition: elementpointlistbase.hh:242
const LocalCoordinateType & localPoint(size_t i) const
obtain local coordinates of i-th integration point
Definition: elementpointlistbase.hh:267
const IntegrationPointListType & quadImp() const
obtain the actual implementation of the quadrature
Definition: elementpointlistbase.hh:342
ElementPointListBase(const GeometryType &elementGeo, const GeometryType &faceGeo, const int localFaceIndex, const QuadratureKeyType &quadKey)
constructor
Definition: elementpointlistbase.hh:226
static const int dimension
dimension of the grid
Definition: elementpointlistbase.hh:210
size_t nop() const
obtain the number of integration points
Definition: elementpointlistbase.hh:252
IntegrationTraits::IntegrationPointListType IntegrationPointListType
type of the integration point list
Definition: elementpointlistbase.hh:213
GeometryType geometry() const
obtain GeometryType for this integration point list
Definition: elementpointlistbase.hh:288
static const int codimension
codimension of the element integration point list
Definition: elementpointlistbase.hh:204
size_t id() const
obtain the identifier of the integration point list
Definition: elementpointlistbase.hh:274
GridPartType::ctype RealType
coordinate type
Definition: elementpointlistbase.hh:207
int order() const
obtain order of the integration point list
Definition: elementpointlistbase.hh:281
GeometryType elementGeometry() const
obtain GeometryType of the corresponding codim-0 the integration point list belongs to
Definition: elementpointlistbase.hh:306
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
constexpr unsigned int dim() const
Return dimension of the type.
Definition: type.hh:360
constexpr unsigned int id() const
Return the topology id of the type.
Definition: type.hh:365
constexpr bool isNone() const
Return true if entity is a singular of any dimension.
Definition: type.hh:355
GeometryType geometryType(const Dune::GeometryType &t)
mapping from GeometryType to VTKGeometryType
Definition: common.hh:151
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
constexpr GeometryType cube(unsigned int dim)
Returns a GeometryType representing a hypercube of dimension dim.
Definition: type.hh:462
constexpr GeometryType none(unsigned int dim)
Returns a GeometryType representing a singular of dimension dim.
Definition: type.hh:471
constexpr GeometryType simplex(unsigned int dim)
Returns a GeometryType representing a simplex of dimension dim.
Definition: type.hh:453
Dune namespace.
Definition: alignedallocator.hh:13
specialize with 'true' for if the codimension 0 entity of the grid part has only one possible geometr...
Definition: capabilities.hh:29
Class providing access to the singletons of the reference elements.
Definition: referenceelements.hh:128
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)