Dune Core Modules (2.8.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_COMMON_GEOMETRY_HH
4#define DUNE_GRID_COMMON_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 //
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 >
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
104 typedef decltype(std::declval<Implementation>().volume()) Volume;
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 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
283 const Geometry &operator= ( const Geometry &rhs ) = delete;
284
285 protected:
286
287 Implementation realGeometry;
288 };
289
290
291
292 //************************************************************************
293 // GEOMETRY Default Implementations
294 //*************************************************************************
295 //
296 // --GeometryDefault
297 //
299 template<int mydim, int cdim, class GridImp, template<int,int,class> class GeometryImp>
301 {
302 public:
303 static const int mydimension = mydim;
304 static const int coorddimension = cdim;
305
306 // save typing
307 typedef typename GridImp::ctype ctype;
308
311
313 typedef ctype Volume;
314
317
320
322 Volume volume () const
323 {
324 GeometryType type = asImp().type();
325
326 // get corresponding reference element
327 auto refElement = referenceElement< ctype , mydim >(type);
328
329 LocalCoordinate localBaryCenter ( 0 );
330 // calculate local bary center
331 const int corners = refElement.size(0,0,mydim);
332 for(int i=0; i<corners; ++i) localBaryCenter += refElement.position(i,mydim);
333 localBaryCenter *= (ctype) (1.0/corners);
334
335 // volume is volume of reference element times integrationElement
336 return refElement.volume() * asImp().integrationElement(localBaryCenter);
337 }
338
341 {
342 GeometryType type = asImp().type();
343
344 // get corresponding reference element
345 auto refElement = referenceElement< ctype , mydim >(type);
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 typedef ctype Volume;
374
376 typedef FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed;
377
379 typedef FieldMatrix< ctype, mydim, cdim > JacobianTransposed;
380
382 FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
383 {
384 return asImp().corner(0);
385 }
386
388 FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& ) const
389 {
390 return FieldVector<ctype, mydim>();
391 }
392
394 Volume volume () const
395 {
396 return Volume(1.0);
397 }
398
400 FieldVector<ctype, cdim> center () const
401 {
402 return asImp().corner(0);
403 }
404
405 private:
406 // Barton-Nackman trick
407 GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
408 const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
409 }; // end GeometryDefault
410
411
412
413 // referenceElement
414 // ----------------
415
416 template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
417 auto referenceElement(const Geometry<mydim,cdim,GridImp,GeometryImp>& geo)
418 -> decltype(referenceElement(geo,geo.impl()))
419 {
420 return referenceElement(geo,geo.impl());
421 }
422
424
439 template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp, typename Impl>
441 -> decltype(referenceElement<typename GridImp::ctype,mydim>(geo.type()))
442 {
444 return referenceElement<typename Geo::ctype,Geo::mydimension>(geo.type());
445 }
446
447} // namespace Dune
448
449#endif // DUNE_GRID_COMMON_GEOMETRY_HH
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
Default implementation for class Geometry.
Definition: geometry.hh:301
Volume volume() const
return volume of the geometry
Definition: geometry.hh:322
GlobalCoordinate center() const
return center of the geometry
Definition: geometry.hh:340
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:319
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:316
ctype Volume
Number type used for the geometry volume.
Definition: geometry.hh:313
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:123
Wrapper class for geometries.
Definition: geometry.hh:67
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:239
Implementation & impl()
access to the underlying implementation
Definition: geometry.hh:81
GeometryType type() const
Return the type 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
@ coorddimension
Definition: geometry.hh:92
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
Volume volume() const
return volume of geometry
Definition: geometry.hh:208
decltype(std::declval< Implementation >().volume()) Volume
Number type used for the geometry volume.
Definition: geometry.hh:104
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:202
const Implementation & impl() const
access to the underlying implementation
Definition: geometry.hh:87
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: geometry.hh:98
GlobalCoordinate center() const
return center of geometry
Definition: geometry.hh:223
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: geometry.hh:134
auto referenceElement(const Geometry< mydim, cdim, GridImp, GeometryImp > &geo, const 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:440
@ mydimension
Definition: geometry.hh:90
const Geometry & operator=(const Geometry &rhs)=delete
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:265
Traits for type conversions and type information.
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:11
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)