Dune Core Modules (2.6.0)

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>
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 
81  Implementation &impl () { return realGeometry; }
87  const Implementation &impl () const { return realGeometry; }
88 
90  enum { mydimension=mydim };
92  enum { coorddimension=cdim };
93 
95  typedef typename GridImp::ctype ctype;
96 
99 
102 
112  typedef typename Implementation::JacobianInverseTransposed JacobianInverseTransposed;
113 
123  typedef typename Implementation::JacobianTransposed JacobianTransposed;
124 
128  GeometryType type () const { return impl().type(); }
129 
131  bool affine() const { return impl().affine(); }
132 
139  int corners () const { return impl().corners(); }
140 
153  GlobalCoordinate corner ( int i ) const
154  {
155  return impl().corner( i );
156  }
157 
163  {
164  return impl().global( local );
165  }
166 
172  {
173  return impl().local( global );
174  }
175 
200  {
201  return impl().integrationElement( local );
202  }
203 
205  ctype volume () const
206  {
207  return impl().volume();
208  }
209 
221  {
222  return impl().center();
223  }
224 
237  {
238  return impl().jacobianTransposed( local );
239  }
240 
263  {
264  return impl().jacobianInverseTransposed(local);
265  }
266  //===========================================================
270  //===========================================================
271 
273  explicit Geometry ( const Implementation &impl )
274  : realGeometry( impl )
275  {}
276 
278 
279  private:
281  const Geometry &operator= ( const Geometry &rhs );
282 
283  protected:
284 
285  Implementation realGeometry;
286  };
287 
288 
289 
290  //************************************************************************
291  // GEOMETRY Default Implementations
292  //*************************************************************************
293  //
294  // --GeometryDefault
295  //
297  template<int mydim, int cdim, class GridImp, template<int,int,class> class GeometryImp>
299  {
300  public:
301  static const int mydimension = mydim;
302  static const int coorddimension = cdim;
303 
304  // save typing
305  typedef typename GridImp::ctype ctype;
306 
309 
312 
315 
317  ctype volume () const
318  {
319  GeometryType type = asImp().type();
320 
321  // get corresponding reference element
322  auto refElement = referenceElement< ctype , mydim >(type);
323 
324  LocalCoordinate localBaryCenter ( 0 );
325  // calculate local bary center
326  const int corners = refElement.size(0,0,mydim);
327  for(int i=0; i<corners; ++i) localBaryCenter += refElement.position(i,mydim);
328  localBaryCenter *= (ctype) (1.0/corners);
329 
330  // volume is volume of reference element times integrationElement
331  return refElement.volume() * asImp().integrationElement(localBaryCenter);
332  }
333 
336  {
337  GeometryType type = asImp().type();
338 
339  // get corresponding reference element
340  auto refElement = referenceElement< ctype , mydim >(type);
341 
342  // center is (for now) the centroid of the reference element mapped to
343  // this geometry.
344  return asImp().global(refElement.position(0,0));
345  }
346 
347  private:
349  GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
350  const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
351  }; // end GeometryDefault
352 
353  template<int cdim, class GridImp, template<int,int,class> class GeometryImp>
354  class GeometryDefaultImplementation<0,cdim,GridImp,GeometryImp>
355  {
356  // my dimension
357  enum { mydim = 0 };
358 
359  public:
360  static const int mydimension = mydim;
361  static const int coorddimension = cdim;
362 
363  // save typing
364  typedef typename GridImp::ctype ctype;
365 
366  typedef FieldVector< ctype, mydim > LocalCoordinate;
367  typedef FieldVector< ctype, cdim > GlobalCoordinate;
368 
370  typedef FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed;
371 
373  typedef FieldMatrix< ctype, mydim, cdim > JacobianTransposed;
374 
376  FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
377  {
378  return asImp().corner(0);
379  }
380 
382  FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& ) const
383  {
384  return FieldVector<ctype, mydim>();
385  }
386 
388  ctype volume () const
389  {
390  return 1.0;
391  }
392 
394  FieldVector<ctype, cdim> center () const
395  {
396  return asImp().corner(0);
397  }
398 
399  private:
400  // Barton-Nackman trick
401  GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
402  const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
403  }; // end GeometryDefault
404 
405 
406 
407  // referenceElement
408  // ----------------
409 
410  template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
411  auto referenceElement(const Geometry<mydim,cdim,GridImp,GeometryImp>& geo)
412  -> decltype(referenceElement(geo,geo.impl()))
413  {
414  return referenceElement(geo,geo.impl());
415  }
416 
418 
433  template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp, typename Impl>
435  -> decltype(referenceElement<typename GridImp::ctype,mydim>(geo.type()))
436  {
438  return referenceElement<typename Geo::ctype,Geo::mydimension>(geo.type());
439  }
440 
441 } // namespace Dune
442 
443 #endif // DUNE_GRID_GEOMETRY_HH
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Default implementation for class Geometry.
Definition: geometry.hh:299
ctype volume() const
return volume of the geometry
Definition: geometry.hh:317
GlobalCoordinate center() const
return center of the geometry
Definition: geometry.hh:335
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:314
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:311
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:277
Wrapper class for geometries.
Definition: geometry.hh:67
const Implementation & impl() const
access to the underlying implementation
Definition: geometry.hh:87
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: geometry.hh:101
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Return the transposed of the Jacobian.
Definition: geometry.hh:236
GeometryType type() const
Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: geometry.hh:128
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the map .
Definition: geometry.hh:162
Geometry(const Implementation &impl)
copy constructor from implementation
Definition: geometry.hh:273
ctype volume() const
return volume of geometry
Definition: geometry.hh:205
Implementation::JacobianTransposed JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:123
@ mydimension
Definition: geometry.hh:90
Implementation::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:112
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: geometry.hh:153
int corners() const
Return the number of corners of the reference element.
Definition: geometry.hh:139
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:95
ctype integrationElement(const LocalCoordinate &local) const
Return the factor appearing in the integral transformation formula.
Definition: geometry.hh:199
Implementation & impl()
access to the underlying implementation
Definition: geometry.hh:81
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: geometry.hh:98
GlobalCoordinate center() const
return center of geometry
Definition: geometry.hh:220
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: geometry.hh:131
@ coorddimension
Definition: geometry.hh:92
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: geometry.hh:434
GeometryImp< mydim, cdim, GridImp > Implementation
type of underlying implementation
Definition: geometry.hh:74
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Return inverse of transposed of Jacobian.
Definition: geometry.hh:262
Implements a matrix constructed from a given type representing a field and compile-time given number ...
unspecified value type referenceElement(T &&... t)
Returns a reference element for the objects t....
Dune namespace.
Definition: alignedallocator.hh:10
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Apr 18, 22:30, 2024)