dune-grid  2.7.1
common/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 
12 #include <dune/common/fmatrix.hh>
13 #include <dune/common/typetraits.hh>
14 
15 #include <dune/geometry/referenceelements.hh>
16 
17 namespace 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  //
30  // Geometry
31  // forwards the interface to the implementation
32  //
33  //*****************************************************************************
34 
65  template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp >
66  class Geometry
67  {
68  public:
74  typedef GeometryImp< mydim, cdim, GridImp > Implementation;
75 
87  const Implementation &impl () const { return realGeometry; }
88 
90  enum { mydimension=mydim };
92  enum { coorddimension=cdim };
93 
95  typedef typename GridImp::ctype ctype;
96 
98  typedef FieldVector<ctype, mydim> LocalCoordinate;
99 
101  typedef FieldVector< ctype, cdim > GlobalCoordinate;
102 
104  typedef decltype(std::declval<Implementation>().volume()) Volume;
105 
116 
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  Volume 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 
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 
310  typedef FieldVector< ctype, mydim > LocalCoordinate;
311  typedef FieldVector< ctype, cdim > GlobalCoordinate;
312 
314  typedef ctype Volume;
315 
317  typedef FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed;
318 
320  typedef FieldMatrix< ctype, mydim, cdim > JacobianTransposed;
321 
323  Volume volume () const
324  {
325  GeometryType type = asImp().type();
326 
327  // get corresponding reference element
328  auto refElement = referenceElement< ctype , mydim >(type);
329 
330  LocalCoordinate localBaryCenter ( 0 );
331  // calculate local bary center
332  const int corners = refElement.size(0,0,mydim);
333  for(int i=0; i<corners; ++i) localBaryCenter += refElement.position(i,mydim);
334  localBaryCenter *= (ctype) (1.0/corners);
335 
336  // volume is volume of reference element times integrationElement
337  return refElement.volume() * asImp().integrationElement(localBaryCenter);
338  }
339 
342  {
343  GeometryType type = asImp().type();
344 
345  // get corresponding reference element
346  auto refElement = referenceElement< ctype , mydim >(type);
347 
348  // center is (for now) the centroid of the reference element mapped to
349  // this geometry.
350  return asImp().global(refElement.position(0,0));
351  }
352 
353  private:
355  GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
356  const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
357  }; // end GeometryDefault
358 
359  template<int cdim, class GridImp, template<int,int,class> class GeometryImp>
360  class GeometryDefaultImplementation<0,cdim,GridImp,GeometryImp>
361  {
362  // my dimension
363  enum { mydim = 0 };
364 
365  public:
366  static const int mydimension = mydim;
367  static const int coorddimension = cdim;
368 
369  // save typing
370  typedef typename GridImp::ctype ctype;
371 
372  typedef FieldVector< ctype, mydim > LocalCoordinate;
373  typedef FieldVector< ctype, cdim > GlobalCoordinate;
374  typedef ctype Volume;
375 
377  typedef FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed;
378 
380  typedef FieldMatrix< ctype, mydim, cdim > JacobianTransposed;
381 
383  FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
384  {
385  return asImp().corner(0);
386  }
387 
389  FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& ) const
390  {
391  return FieldVector<ctype, mydim>();
392  }
393 
395  Volume volume () const
396  {
397  return Volume(1.0);
398  }
399 
401  FieldVector<ctype, cdim> center () const
402  {
403  return asImp().corner(0);
404  }
405 
406  private:
407  // Barton-Nackman trick
408  GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
409  const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
410  }; // end GeometryDefault
411 
412 
413 
414  // referenceElement
415  // ----------------
416 
417  template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
419  -> decltype(referenceElement(geo,geo.impl()))
420  {
421  return referenceElement(geo,geo.impl());
422  }
423 
425 
440  template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp, typename Impl>
442  -> decltype(referenceElement<typename GridImp::ctype,mydim>(geo.type()))
443  {
445  return referenceElement<typename Geo::ctype,Geo::mydimension>(geo.type());
446  }
447 
448 } // namespace Dune
449 
450 #endif // DUNE_GRID_GEOMETRY_HH
Include standard header files.
Definition: agrid.hh:59
auto referenceElement(const Geometry< mydim, cdim, GridImp, GeometryImp > &geo) -> decltype(referenceElement(geo, geo.impl()))
Definition: common/geometry.hh:418
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:180
Wrapper class for geometries.
Definition: common/geometry.hh:67
const Implementation & impl() const
access to the underlying implementation
Definition: common/geometry.hh:87
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: common/geometry.hh:101
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Return the transposed of the Jacobian.
Definition: common/geometry.hh:239
GeometryType type() const
Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: common/geometry.hh:131
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the map .
Definition: common/geometry.hh:165
Geometry(const Implementation &impl)
copy constructor from implementation
Definition: common/geometry.hh:276
@ mydimension
Definition: common/geometry.hh:90
Implementation::JacobianTransposed JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:126
Implementation::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:115
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: common/geometry.hh:156
decltype(std::declval< Implementation >().volume()) typedef Volume
Number type used for the geometry volume.
Definition: common/geometry.hh:104
int corners() const
Return the number of corners of the reference element.
Definition: common/geometry.hh:142
Volume volume() const
return volume of geometry
Definition: common/geometry.hh:208
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: common/geometry.hh:95
ctype integrationElement(const LocalCoordinate &local) const
Return the factor appearing in the integral transformation formula.
Definition: common/geometry.hh:202
Implementation & impl()
access to the underlying implementation
Definition: common/geometry.hh:81
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: common/geometry.hh:98
Implementation realGeometry
Definition: common/geometry.hh:288
GlobalCoordinate center() const
return center of geometry
Definition: common/geometry.hh:223
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: common/geometry.hh:134
auto referenceElement(const Geometry< mydim, cdim, GridImp, GeometryImp > &geo, const Impl &impl) -> decltype(referenceElement< typename GridImp::ctype, mydim >(geo.type()))
Second-level dispatch to select the correct reference element for a grid geometry.
Definition: common/geometry.hh:441
@ coorddimension
Definition: common/geometry.hh:92
GeometryImp< mydim, cdim, GridImp > Implementation
type of underlying implementation
Definition: common/geometry.hh:74
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Return inverse of transposed of Jacobian.
Definition: common/geometry.hh:265
Default implementation for class Geometry.
Definition: common/geometry.hh:302
FieldVector< ctype, cdim > GlobalCoordinate
Definition: common/geometry.hh:311
static const int mydimension
Definition: common/geometry.hh:304
Volume volume() const
return volume of the geometry
Definition: common/geometry.hh:323
FieldVector< ctype, mydim > LocalCoordinate
Definition: common/geometry.hh:310
GlobalCoordinate center() const
return center of the geometry
Definition: common/geometry.hh:341
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:320
GridImp::ctype ctype
Definition: common/geometry.hh:308
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:317
static const int coorddimension
Definition: common/geometry.hh:305
ctype Volume
Number type used for the geometry volume.
Definition: common/geometry.hh:314
GridImp::ctype ctype
Definition: common/geometry.hh:370
FieldVector< ctype, cdim > GlobalCoordinate
Definition: common/geometry.hh:373
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:377
FieldVector< ctype, cdim > global(const FieldVector< ctype, mydim > &local) const
return the only coordinate
Definition: common/geometry.hh:383
FieldVector< ctype, mydim > LocalCoordinate
Definition: common/geometry.hh:372
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:380
FieldVector< ctype, cdim > center() const
return center of the geometry
Definition: common/geometry.hh:401
Volume volume() const
return volume of the geometry
Definition: common/geometry.hh:395