Dune Core Modules (2.3.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 namespace FacadeOptions
29 {
30
33
45 template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp >
47 {
49 static const bool v = true;
50 };
51
52 }
53
54
55
56 //*****************************************************************************
57 //
58 // Geometry
59 // forwards the interface to the implementation
60 //
61 //*****************************************************************************
62
100 template< int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp >
102 {
103 #if DUNE_GRID_EXPERIMENTAL_GRID_EXTENSIONS
104 public:
105 #else
106 protected:
107 // give the GridDefaultImplementation class access to the realImp
108 friend class GridDefaultImplementation<
109 GridImp::dimension, GridImp::dimensionworld,
110 typename GridImp::ctype,
111 typename GridImp::GridFamily> ;
112 #endif
113 // type of underlying implementation, for internal use only
114 typedef GeometryImp< mydim, cdim, GridImp > Implementation;
115
117 const Implementation &impl () const { return realGeometry; }
118
119 public:
121 enum { dimension=GridImp::dimension };
123 enum { mydimension=mydim };
125 enum { coorddimension=cdim };
126
128 enum { dimensionworld=GridImp::dimensionworld };
130 typedef typename GridImp::ctype ctype;
131
134
137
139 typedef typename Implementation::JacobianInverseTransposed JacobianInverseTransposed;
140
143 typedef JacobianInverseTransposed Jacobian DUNE_DEPRECATED_MSG ( "type Geometry::Jacobian is deprecated, use Geometry::JacobianInverseTransposed instead." );
144
146 typedef typename Implementation::JacobianTransposed JacobianTransposed;
147
151 GeometryType type () const { return impl().type(); }
152
154 bool affine() const { return impl().affine(); }
155
162 int corners () const { return impl().corners(); }
163
176 GlobalCoordinate corner ( int i ) const
177 {
178 return impl().corner( i );
179 }
180
186 {
187 return impl().global( local );
188 }
189
195 {
196 return impl().local( global );
197 }
198
223 {
224 return impl().integrationElement( local );
225 }
226
228 ctype volume () const
229 {
230 return impl().volume();
231 }
232
244 {
245 return impl().center();
246 }
247
257 const JacobianTransposed &
258 jacobianTransposed ( const LocalCoordinate& local ) const
259 {
260 return impl().jacobianTransposed( local );
261 }
262
284 {
285 return impl().jacobianInverseTransposed(local);
286 }
287
289 explicit Geometry ( const Implementation &impl )
290 : realGeometry( impl )
291 {
292 deprecationWarning ( integral_constant< bool, storeReference >() );
293 }
294
295 private:
297 const Geometry &operator= ( const Geometry &rhs );
298
299 void DUNE_DEPRECATED_MSG( "This Dune::Geometry is still a reference to its implementation." )
300 deprecationWarning ( integral_constant< bool, true > ) {}
301
302 void
303 deprecationWarning ( integral_constant< bool, false > ) {}
304
305 protected:
307
309 };
310
311
312
313 //************************************************************************
314 // GEOMETRY Default Implementations
315 //*************************************************************************
316 //
317 // --GeometryDefault
318 //
320 template<int mydim, int cdim, class GridImp, template<int,int,class> class GeometryImp>
322 {
323 public:
324 static const int mydimension = mydim;
325 static const int coorddimension = cdim;
326
327 // save typing
328 typedef typename GridImp::ctype ctype;
329
332
335
338
340 ctype volume () const
341 {
342 GeometryType type = asImp().type();
343
344 // get corresponding reference element
345 const ReferenceElement< ctype , mydim > & refElement =
347
348 LocalCoordinate localBaryCenter ( 0 );
349 // calculate local bary center
350 const int corners = refElement.size(0,0,mydim);
351 for(int i=0; i<corners; ++i) localBaryCenter += refElement.position(i,mydim);
352 localBaryCenter *= (ctype) (1.0/corners);
353
354 // volume is volume of reference element times integrationElement
355 return refElement.volume() * asImp().integrationElement(localBaryCenter);
356 }
357
360 {
361 GeometryType type = asImp().type();
362
363 // get corresponding reference element
364 const ReferenceElement< ctype , mydim > & refElement =
366
367 // center is (for now) the centroid of the reference element mapped to
368 // this geometry.
369 return asImp().global(refElement.position(0,0));
370 }
371
372 private:
374 GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
375 const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
376 }; // end GeometryDefault
377
378 template<int cdim, class GridImp, template<int,int,class> class GeometryImp>
379 class GeometryDefaultImplementation<0,cdim,GridImp,GeometryImp>
380 {
381 // my dimension
382 enum { mydim = 0 };
383
384 public:
385 static const int mydimension = mydim;
386 static const int coorddimension = cdim;
387
388 // save typing
389 typedef typename GridImp::ctype ctype;
390
391 typedef FieldVector< ctype, mydim > LocalCoordinate;
392 typedef FieldVector< ctype, cdim > GlobalCoordinate;
393
396
399
401 FieldVector<ctype, cdim> global (const FieldVector<ctype, mydim>& local) const
402 {
403 return asImp().corner(0);
404 }
405
407 FieldVector<ctype, mydim> local (const FieldVector<ctype, cdim>& ) const
408 {
409 return FieldVector<ctype, mydim>();
410 }
411
413 bool checkInside (const FieldVector<ctype, mydim>& ) const
414 {
415 return true;
416 }
417
419 ctype volume () const
420 {
421 return 1.0;
422 }
423
425 FieldVector<ctype, cdim> center () const
426 {
427 return asImp().corner(0);
428 }
429
430 private:
431 // Barton-Nackman trick
432 GeometryImp<mydim,cdim,GridImp>& asImp () {return static_cast<GeometryImp<mydim,cdim,GridImp>&>(*this);}
433 const GeometryImp<mydim,cdim,GridImp>& asImp () const {return static_cast<const GeometryImp<mydim,cdim,GridImp>&>(*this);}
434 }; // end GeometryDefault
435
436} // namespace Dune
437
438#endif // DUNE_GRID_GEOMETRY_HH
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Default implementation for class Geometry.
Definition: geometry.hh:322
ctype volume() const
return volume of the geometry
Definition: geometry.hh:340
GlobalCoordinate center() const
return center of the geometry
Definition: geometry.hh:359
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:337
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:334
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
Wrapper class for geometries.
Definition: geometry.hh:102
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: geometry.hh:136
@ dimension
Definition: geometry.hh:121
GeometryType type() const
Return the name of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: geometry.hh:151
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the map .
Definition: geometry.hh:185
@ dimensionworld
Definition: geometry.hh:128
@ coorddimension
Definition: geometry.hh:125
Geometry(const Implementation &impl)
copy constructor from implementation
Definition: geometry.hh:289
ctype volume() const
return volume of geometry
Definition: geometry.hh:228
Implementation::JacobianTransposed JacobianTransposed
type of jacobian transposed
Definition: geometry.hh:146
const JacobianInverseTransposed & jacobianInverseTransposed(const LocalCoordinate &local) const
Return inverse of transposed of Jacobian.
Definition: geometry.hh:283
Implementation::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian inverse transposed
Definition: geometry.hh:139
const JacobianTransposed & jacobianTransposed(const LocalCoordinate &local) const
Return the transposed of the Jacobian.
Definition: geometry.hh:258
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: geometry.hh:176
int corners() const
Return the number of corners of the reference element.
Definition: geometry.hh:162
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: geometry.hh:130
ctype integrationElement(const LocalCoordinate &local) const
Return the factor appearing in the integral transformation formula.
Definition: geometry.hh:222
const Implementation & impl() const
return reference to the implementation
Definition: geometry.hh:117
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: geometry.hh:133
GlobalCoordinate center() const
return center of geometry
Definition: geometry.hh:243
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: geometry.hh:154
JacobianInverseTransposed Jacobian DUNE_DEPRECATED_MSG("type Geometry::Jacobian is deprecated, use Geometry::JacobianInverseTransposed instead.")
Definition: geometry.hh:143
@ mydimension
Definition: geometry.hh:123
Definition: grid.hh:1017
This class provides access to geometric and topological properties of a reference element....
Definition: referenceelements.hh:58
int size(int c) const
number of subentities of codimension c
Definition: referenceelements.hh:86
const FieldVector< ctype, dim > & position(int i, int c) const
position of the barycenter of entity (i,c)
Definition: referenceelements.hh:154
ctype volume() const
obtain the volume of the reference element
Definition: referenceelements.hh:280
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Dune namespace.
Definition: alignment.hh:14
Traits class determining whether the Dune::Geometry facade class stores the implementation object by ...
Definition: geometry.hh:47
static const bool v
Whether to store by reference.
Definition: geometry.hh:49
static const ReferenceElement< ctype, dim > & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:568
T1 type
The selected type.
Definition: typetraits.hh:426
Generate a type for a given integral constant.
Definition: typetraits.hh:457
Traits for type conversions and type information.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)