Dune Core Modules (2.5.0)

geometry.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_ALBERTA_GEOMETRY_HH
4#define DUNE_ALBERTA_GEOMETRY_HH
5
7#include <dune/grid/albertagrid/misc.hh>
9
10#if HAVE_ALBERTA
11
12namespace Dune
13{
14
15 // Forward Declarations
16 // --------------------
17
18 template< int dim, int dimworld >
19 class AlbertaGrid;
20
21
22
23 // AlbertaGridCoordinateReader
24 // ---------------------------
25
26 template< int codim, class GridImp >
27 struct AlbertaGridCoordinateReader
28 {
29 typedef typename std::remove_const< GridImp >::type Grid;
30
31 static const int dimension = Grid::dimension;
32 static const int codimension = codim;
33 static const int mydimension = dimension - codimension;
34 static const int coorddimension = Grid::dimensionworld;
35
36 typedef Alberta::Real ctype;
37
38 typedef Alberta::ElementInfo< dimension > ElementInfo;
39 typedef FieldVector< ctype, coorddimension > Coordinate;
40
41 AlbertaGridCoordinateReader ( const GridImp &grid,
42 const ElementInfo &elementInfo,
43 int subEntity )
44 : grid_( grid ),
45 elementInfo_( elementInfo ),
46 subEntity_( subEntity )
47 {}
48
49 const ElementInfo &elementInfo () const
50 {
51 return elementInfo_;
52 }
53
54 void coordinate ( int i, Coordinate &x ) const
55 {
56 assert( !elementInfo_ == false );
57 assert( (i >= 0) && (i <= mydimension) );
58
59 const int k = mapVertices( subEntity_, i );
60 const Alberta::GlobalVector &coord = grid_.getCoord( elementInfo_, k );
61 for( int j = 0; j < coorddimension; ++j )
62 x[ j ] = coord[ j ];
63 }
64
65 bool hasDeterminant () const
66 {
67 return false;
68 }
69
70 ctype determinant () const
71 {
72 assert( hasDeterminant() );
73 return ctype( 0 );
74 }
75
76 private:
77 static int mapVertices ( int subEntity, int i )
78 {
79 return Alberta::MapVertices< dimension, codimension >::apply( subEntity, i );
80 }
81
82 const Grid &grid_;
83 const ElementInfo &elementInfo_;
84 const int subEntity_;
85 };
86
87
88
89 // AlbertaGridGeometry
90 // -------------------
91
104 template< int mydim, int cdim, class GridImp >
106 {
108
109 // remember type of the grid
110 typedef GridImp Grid;
111
112 // dimension of barycentric coordinates
113 static const int dimbary = mydim + 1;
114
115 public:
117 typedef Alberta::Real ctype;
118
119 static const int dimension = Grid :: dimension;
120 static const int mydimension = mydim;
121 static const int codimension = dimension - mydimension;
122 static const int coorddimension = cdim;
123
126
129
130 private:
131 static const int numCorners = mydimension + 1;
132
134
135 public:
137 {
138 invalidate();
139 }
140
141 template< class CoordReader >
142 AlbertaGridGeometry ( const CoordReader &coordReader )
143 {
144 build( coordReader );
145 }
146
149 {
150 typedef typename Impl::SimplexTopology< mydimension >::type Topology;
151 return GeometryType( Topology() );
152 }
153
155 bool affine () const { return true; }
156
158 int corners () const
159 {
160 return numCorners;
161 }
162
164 GlobalCoordinate corner ( const int i ) const
165 {
166 assert( (i >= 0) && (i < corners()) );
167 return coord_[ i ];
168 }
169
172 {
173 return centroid_;
174 }
175
177 GlobalCoordinate global ( const LocalCoordinate &local ) const;
178
181
188 {
189 assert( calcedDet_ );
190 return elDet_;
191 }
192
195 {
196 return integrationElement();
197 }
198
200 ctype volume () const
201 {
203 }
204
210 const JacobianTransposed &jacobianTransposed () const;
211
213 const JacobianTransposed &
214 jacobianTransposed ( const LocalCoordinate &local ) const
215 {
216 return jacobianTransposed();
217 }
218
224 const JacobianInverseTransposed &jacobianInverseTransposed () const;
225
227 const JacobianInverseTransposed &
229 {
231 }
232
233 //***********************************************************************
234 // Methods that not belong to the Interface, but have to be public
235 //***********************************************************************
236
239 {
240 builtJT_ = false;
241 builtJTInv_ = false;
242 calcedDet_ = false;
243 }
244
246 template< class CoordReader >
247 void build ( const CoordReader &coordReader );
248
249 void print ( std::ostream &out ) const;
250
251 private:
252 // calculates the volume of the element
253 ctype elDeterminant () const
254 {
255 return std::abs( Alberta::determinant( jacobianTransposed() ) );
256 }
257
259 CoordMatrix coord_;
260
262 GlobalCoordinate centroid_;
263
264 // storage for the transposed of the jacobian
265 mutable JacobianTransposed jT_;
266
267 // storage for the transposed inverse of the jacboian
268 mutable JacobianInverseTransposed jTInv_;
269
270 // has jT_ been computed, yet?
271 mutable bool builtJT_;
272 // has jTInv_ been computed, yet?
273 mutable bool builtJTInv_;
274
275 mutable bool calcedDet_;
276 mutable ctype elDet_;
277 };
278
279
280
281 // AlbertaGridGlobalGeometry
282 // -------------------------
283
284 template< int mydim, int cdim, class GridImp >
285 class AlbertaGridGlobalGeometry
286 : public AlbertaGridGeometry< mydim, cdim, GridImp >
287 {
288 typedef AlbertaGridGlobalGeometry< mydim, cdim, GridImp > This;
289 typedef AlbertaGridGeometry< mydim, cdim, GridImp > Base;
290
291 public:
292 AlbertaGridGlobalGeometry ()
293 : Base()
294 {}
295
296 template< class CoordReader >
297 AlbertaGridGlobalGeometry ( const CoordReader &coordReader )
298 : Base( coordReader )
299 {}
300 };
301
302
303#if !DUNE_ALBERTA_CACHE_COORDINATES
304 template< int dim, int cdim >
305 class AlbertaGridGlobalGeometry< dim, cdim, const AlbertaGrid< dim, cdim > >
306 {
307 typedef AlbertaGridGlobalGeometry< dim, cdim, const AlbertaGrid< dim, cdim > > This;
308
309 // remember type of the grid
310 typedef AlbertaGrid< dim, cdim > Grid;
311
312 // dimension of barycentric coordinates
313 static const int dimbary = dim + 1;
314
315 typedef Alberta::ElementInfo< dim > ElementInfo;
316
317 public:
319 typedef Alberta::Real ctype;
320
321 static const int dimension = Grid::dimension;
322 static const int mydimension = dimension;
323 static const int codimension = dimension - mydimension;
324 static const int coorddimension = cdim;
325
326 typedef FieldVector< ctype, mydimension > LocalCoordinate;
327 typedef FieldVector< ctype, coorddimension > GlobalCoordinate;
328
329 typedef FieldMatrix< ctype, mydimension, coorddimension > JacobianTransposed;
330 typedef FieldMatrix< ctype, coorddimension, mydimension > JacobianInverseTransposed;
331
332 private:
333 static const int numCorners = mydimension + 1;
334
336
337 public:
338 AlbertaGridGlobalGeometry ()
339 {
340 invalidate();
341 }
342
343 template< class CoordReader >
344 AlbertaGridGlobalGeometry ( const CoordReader &coordReader )
345 {
346 build( coordReader );
347 }
348
350 GeometryType type () const
351 {
352 typedef typename Impl::SimplexTopology< mydimension >::type Topology;
353 return GeometryType( Topology() );
354 }
355
357 int corners () const
358 {
359 return numCorners;
360 }
361
363 GlobalCoordinate corner ( const int i ) const
364 {
365 assert( (i >= 0) && (i < corners()) );
366 const Alberta::GlobalCoordinate &x = elementInfo_.coordinate( i );
367 GlobalCoordinate y;
368 for( int j = 0; j < coorddimension; ++j )
369 y[ j ] = x[ j ];
370 return y;
371 }
372
374 GlobalCoordinate center () const
375 {
376 GlobalCoordinate centroid_ = corner( 0 );
377 for( int i = 1; i < numCorners; ++i )
378 centroid_ += corner( i );
379 centroid_ *= ctype( 1 ) / ctype( numCorners );
380 return centroid_;
381 }
382
384 GlobalCoordinate global ( const LocalCoordinate &local ) const;
385
387 LocalCoordinate local ( const GlobalCoordinate &global ) const;
388
394 ctype integrationElement () const
395 {
396 return elementInfo_.geometryCache().integrationElement();
397 }
398
400 ctype integrationElement ( const LocalCoordinate &local ) const
401 {
402 return integrationElement();
403 }
404
406 ctype volume () const
407 {
408 return integrationElement() / ctype( Factorial< mydimension >::factorial );
409 }
410
416 const JacobianTransposed &jacobianTransposed () const
417 {
418 return elementInfo_.geometryCache().jacobianTransposed();
419 }
420
422 const JacobianTransposed &
423 jacobianTransposed ( const LocalCoordinate &local ) const
424 {
425 return jacobianTransposed();
426 }
427
433 const JacobianInverseTransposed &jacobianInverseTransposed () const
434 {
435 return elementInfo_.geometryCache().jacobianInverseTransposed();
436 }
437
439 const JacobianInverseTransposed &
440 jacobianInverseTransposed ( const LocalCoordinate &local ) const
441 {
443 }
444
445 //***********************************************************************
446 // Methods that not belong to the Interface, but have to be public
447 //***********************************************************************
448
450 void invalidate ()
451 {
452 elementInfo_ = ElementInfo();
453 }
454
456 template< class CoordReader >
457 void build ( const CoordReader &coordReader )
458 {
459 elementInfo_ = coordReader.elementInfo();
460 }
461
462 private:
463 ElementInfo elementInfo_;
464 };
465#endif // #if !DUNE_ALBERTA_CACHE_COORDINATES
466
467
468
469 // AlbertaGridLocalGeometryProvider
470 // --------------------------------
471
472 template< class Grid >
473 class AlbertaGridLocalGeometryProvider
474 {
475 typedef AlbertaGridLocalGeometryProvider< Grid > This;
476
477 public:
478 typedef typename Grid::ctype ctype;
479
480 static const int dimension = Grid::dimension;
481
482 template< int codim >
483 struct Codim
484 {
485 typedef AlbertaGridGeometry< dimension-codim, dimension, Grid > LocalGeometry;
486 };
487
488 typedef typename Codim< 0 >::LocalGeometry LocalElementGeometry;
489 typedef typename Codim< 1 >::LocalGeometry LocalFaceGeometry;
490
491 static const int numChildren = 2;
492 static const int numFaces = dimension + 1;
493
494 static const int minFaceTwist = Alberta::Twist< dimension, dimension-1 >::minTwist;
495 static const int maxFaceTwist = Alberta::Twist< dimension, dimension-1 >::maxTwist;
496 static const int numFaceTwists = maxFaceTwist - minFaceTwist + 1;
497
498 private:
499 struct GeoInFatherCoordReader;
500 struct FaceCoordReader;
501
502 AlbertaGridLocalGeometryProvider ()
503 {
504 buildGeometryInFather();
505 buildFaceGeometry();
506 }
507
508 ~AlbertaGridLocalGeometryProvider ()
509 {
510 for( int child = 0; child < numChildren; ++child )
511 {
512 delete geometryInFather_[ child ][ 0 ];
513 delete geometryInFather_[ child ][ 1 ];
514 }
515
516 for( int i = 0; i < numFaces; ++i )
517 {
518 for( int j = 0; j < numFaceTwists; ++j )
519 delete faceGeometry_[ i ][ j ];
520 }
521 }
522
523 void buildGeometryInFather();
524 void buildFaceGeometry();
525
526 public:
527 const LocalElementGeometry &
528 geometryInFather ( int child, const int orientation = 1 ) const
529 {
530 assert( (child >= 0) && (child < numChildren) );
531 assert( (orientation == 1) || (orientation == -1) );
532 return *geometryInFather_[ child ][ (orientation + 1) / 2 ];
533 }
534
535 const LocalFaceGeometry &
536 faceGeometry ( int face, int twist = 0 ) const
537 {
538 assert( (face >= 0) && (face < numFaces) );
539 assert( (twist >= minFaceTwist) && (twist <= maxFaceTwist) );
540 return *faceGeometry_[ face ][ twist - minFaceTwist ];
541 }
542
543 static const This &instance ()
544 {
545 static This theInstance;
546 return theInstance;
547 }
548
549 private:
550 template< int codim >
551 static int mapVertices ( int subEntity, int i )
552 {
553 return Alberta::MapVertices< dimension, codim >::apply( subEntity, i );
554 }
555
556 const LocalElementGeometry *geometryInFather_[ numChildren ][ 2 ];
557 const LocalFaceGeometry *faceGeometry_[ numFaces ][ numFaceTwists ];
558 };
559
560} // namespace Dune
561
562#endif // #if HAVE_ALBERTA
563
564#endif // #ifndef DUNE_ALBERTA_GEOMETRY_HH
geometry implementation for AlbertaGrid
Definition: geometry.hh:106
GeometryType type() const
obtain the type of reference element
Definition: geometry.hh:148
GlobalCoordinate center() const
return center of geometry
Definition: geometry.hh:171
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
transposed of the geometry mapping's Jacobian
Definition: geometry.hh:214
const JacobianTransposed & jacobianTransposed() const
transposed of the geometry mapping's Jacobian
Definition: geometry.cc:53
int corners() const
number of corner the geometry
Definition: geometry.hh:158
const JacobianInverseTransposed & jacobianInverseTransposed() const
transposed inverse of the geometry mapping's Jacobian
Definition: geometry.cc:71
GlobalCoordinate corner(const int i) const
obtain the i-th corner of this geometry
Definition: geometry.hh:164
ctype integrationElement() const
integration element of the geometry mapping
Definition: geometry.hh:187
ctype integrationElement(const LocalCoordinate &local) const
integration element of the geometry mapping
Definition: geometry.hh:194
GlobalCoordinate global(const LocalCoordinate &local) const
map a point from the refence element to the geometry
Definition: geometry.cc:32
Alberta::Real ctype
type of coordinates
Definition: geometry.hh:117
ctype volume() const
volume of geometry
Definition: geometry.hh:200
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
transposed inverse of the geometry mapping's Jacobian
Definition: geometry.hh:228
void invalidate()
invalidate the geometry
Definition: geometry.hh:238
bool affine() const
returns always true since we only have simplices
Definition: geometry.hh:155
void build(const CoordReader &coordReader)
build the geometry from a coordinate reader
Definition: geometry.cc:88
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:268
@ dimensionworld
The dimension of the world the grid lives in.
Definition: grid.hh:393
@ dimension
The dimension of the grid.
Definition: grid.hh:387
ct ctype
Define type used for coordinates in grid module.
Definition: grid.hh:522
Wrapper and interface classes for element geometries.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
provides a wrapper for ALBERTA's el_info structure
Dune namespace.
Definition: alignment.hh:11
Calculates the factorial of m at compile time.
Definition: math.hh:87
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)