Dune Core Modules (2.4.1)

geometry.hh
Go to the documentation of this file.
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_GRID_GEOMETRY_HH
4#define DUNE_GRID_GEOMETRY_HH
5
10#include <cassert>
11
14
15#include <dune/geometry/referenceelements.hh>
16
17namespace Dune
18{
19
20 // External Forward Declarations
21 // -----------------------------
22
23 template< int dim, int dimworld, class ct, class GridFamily >
24 class GridDefaultImplementation;
25
26
27 //*****************************************************************************
28 //
29 // Geometry
30 // forwards the interface to the implementation
31 //
32 //*****************************************************************************
33
64 template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp >
66 {
67 #if DUNE_GRID_EXPERIMENTAL_GRID_EXTENSIONS
68 public:
69 #else
70 protected:
71 // give the GridDefaultImplementation class access to the realImp
72 friend class GridDefaultImplementation<
73 GridImp::dimension, GridImp::dimensionworld,
74 typename GridImp::ctype,
75 typename GridImp::GridFamily> ;
76 #endif
77 // type of underlying implementation, for internal use only
78 typedef GeometryImp< mydim, cdim, GridImp > Implementation;
79
81 const Implementation &impl () const { return realGeometry; }
82
83 public:
87 enum { dimension=GridImp::dimension };
89 enum { mydimension=mydim };
91 enum { coorddimension=cdim };
92
96 enum { dimensionworld=GridImp::dimensionworld };
98 typedef typename GridImp::ctype ctype;
99
102
105
115 typedef typename Implementation::JacobianInverseTransposed JacobianInverseTransposed;
116
126 typedef typename Implementation::JacobianTransposed JacobianTransposed;
127
131 GeometryType type () const { return impl().type(); }
132
134 bool affine() const { return impl().affine(); }
135
142 int corners () const { return impl().corners(); }
143
156 GlobalCoordinate corner ( int i ) const
157 {
158 return impl().corner( i );
159 }
160
166 {
167 return impl().global( local );
168 }
169
175 {
176 return impl().local( global );
177 }
178
203 {
204 return impl().integrationElement( local );
205 }
206
208 ctype volume () const
209 {
210 return impl().volume();
211 }
212
224 {
225 return impl().center();
226 }
227
240 {
241 return impl().jacobianTransposed( local );
242 }
243
266 {
267 return impl().jacobianInverseTransposed(local);
268 }
269 //===========================================================
273 //===========================================================
274
276 explicit Geometry ( const Implementation &impl )
277 : realGeometry( impl )
278 {}
279
281
282 private:
284 const Geometry &operator= ( const Geometry &rhs );
285
286 protected:
287
288 Implementation realGeometry;
289 };
290
291
292
293 //************************************************************************
294 // GEOMETRY Default Implementations
295 //*************************************************************************
296 //
297 // --GeometryDefault
298 //
300 template<int mydim, int cdim, class GridImp, template<int,int,class> class GeometryImp>
302 {
303 public:
304 static const int mydimension = mydim;
305 static const int coorddimension = cdim;
306
307 // save typing
308 typedef typename GridImp::ctype ctype;
309
312
315
318
320 ctype volume () const
321 {
322 GeometryType type = asImp().type();
323
324 // get corresponding reference element
325 const ReferenceElement< ctype , mydim > & refElement =
327
328 LocalCoordinate localBaryCenter ( 0 );
329 // calculate local bary center
330 const int corners = refElement.size(0,0,mydim);
331 for(int i=0; i<corners; ++i) localBaryCenter += refElement.position(i,mydim);
332 localBaryCenter *= (ctype) (1.0/corners);
333
334 // volume is volume of reference element times integrationElement
335 return refElement.volume() * asImp().integrationElement(localBaryCenter);
336 }
337
340 {
341 GeometryType type = asImp().type();
342
343 // get corresponding reference element
344 const ReferenceElement< ctype , mydim > & refElement =
346
347 // center is (for now) the centroid of the reference element mapped to
348 // this geometry.
349 return asImp().global(refElement.position(0,0));
350 }
351
352 private:
354 GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
355 const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
356 }; // end GeometryDefault
357
358 template<int cdim, class GridImp, template<int,int,class> class GeometryImp>
359 class GeometryDefaultImplementation<0,cdim,GridImp,GeometryImp>
360 {
361 // my dimension
362 enum { mydim = 0 };
363
364 public:
365 static const int mydimension = mydim;
366 static const int coorddimension = cdim;
367
368 // save typing
369 typedef typename GridImp::ctype ctype;
370
371 typedef FieldVector< ctype, mydim > LocalCoordinate;
372 typedef FieldVector< ctype, cdim > GlobalCoordinate;
373
376
379
381 FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
382 {
383 return asImp().corner(0);
384 }
385
387 FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& ) const
388 {
389 return FieldVector<ctype, mydim>();
390 }
391
393 ctype volume () const
394 {
395 return 1.0;
396 }
397
399 FieldVector<ctype, cdim> center () const
400 {
401 return asImp().corner(0);
402 }
403
404 private:
405 // Barton-Nackman trick
406 GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
407 const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
408 }; // end GeometryDefault
409
410} // namespace Dune
411
412#endif // DUNE_GRID_GEOMETRY_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:94
Default implementation for class Geometry.
Definition: geometry.hh:302
ctype volume() const
return volume of the geometry
Definition: geometry.hh:320
GlobalCoordinate center() const
return center of the geometry
Definition: geometry.hh:339
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:317
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:314
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
Wrapper class for geometries.
Definition: geometry.hh:66
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: geometry.hh:104
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Return the transposed of the Jacobian.
Definition: geometry.hh:239
@ coorddimension
Definition: geometry.hh:91
GeometryType type() const
Return the name of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: geometry.hh:131
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the map .
Definition: geometry.hh:165
Geometry(const Implementation &impl)
copy constructor from implementation
Definition: geometry.hh:276
@ mydimension
Definition: geometry.hh:89
ctype volume() const
return volume of geometry
Definition: geometry.hh:208
Implementation::JacobianTransposed JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:126
Implementation::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:115
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: geometry.hh:156
int corners() const
Return the number of corners of the reference element.
Definition: geometry.hh:142
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:98
ctype integrationElement(const LocalCoordinate &local) const
Return the factor appearing in the integral transformation formula.
Definition: geometry.hh:202
const Implementation & impl() const
return reference to the implementation
Definition: geometry.hh:81
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: geometry.hh:101
GlobalCoordinate center() const
return center of geometry
Definition: geometry.hh:223
@ dimension
Definition: geometry.hh:87
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: geometry.hh:134
@ dimensionworld
Definition: geometry.hh:96
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Return inverse of transposed of Jacobian.
Definition: geometry.hh:265
Definition: grid.hh:1030
This class provides access to geometric and topological properties of a reference element.
Definition: referenceelements.hh:55
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:82
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelements.hh:150
ctype volume() const
obtain the volume of the reference element
Definition: referenceelements.hh:210
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Dune namespace.
Definition: alignment.hh:10
static const ReferenceElement< ctype, dim > & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:484
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)