3 #ifndef DUNE_GRID_GEOMETRY_HH 4 #define DUNE_GRID_GEOMETRY_HH 12 #include <dune/common/fmatrix.hh> 13 #include <dune/common/typetraits.hh> 15 #include <dune/geometry/referenceelements.hh> 23 template<
int dim,
int dimworld,
class ct,
class Gr
idFamily >
64 template<
int mydim,
int cdim,
class Gr
idImp,
template<
int,
int,
class >
class GeometryImp >
67 #if DUNE_GRID_EXPERIMENTAL_GRID_EXTENSIONS 73 GridImp::dimension, GridImp::dimensionworld,
74 typename GridImp::ctype,
75 typename GridImp::GridFamily> ;
90 typedef typename GridImp::ctype
ctype;
150 return impl().corner( i );
159 return impl().global( local );
168 return impl().local( global );
196 return impl().integrationElement( local );
202 return impl().volume();
217 return impl().center();
233 return impl().jacobianTransposed( local );
259 return impl().jacobianInverseTransposed(local);
292 template<
int mydim,
int cdim,
class Gr
idImp,
template<
int,
int,
class>
class GeometryImp>
300 typedef typename GridImp::ctype
ctype;
317 const ReferenceElement< ctype , mydim > & refElement =
318 ReferenceElements< ctype, mydim >::general(type);
320 LocalCoordinate localBaryCenter ( 0 );
322 const int corners = refElement.size(0,0,mydim);
323 for(
int i=0; i<
corners; ++i) localBaryCenter += refElement.position(i,mydim);
324 localBaryCenter *= (
ctype) (1.0/corners);
327 return refElement.volume() * asImp().integrationElement(localBaryCenter);
336 const ReferenceElement< ctype , mydim > & refElement =
337 ReferenceElements< ctype, mydim >::general(type);
341 return asImp().global(refElement.position(0,0));
346 GeometryImp<mydim,cdim,GridImp>& asImp () {
return static_cast<GeometryImp<mydim,cdim,GridImp>&
>(*this);}
347 const GeometryImp<mydim,cdim,GridImp>& asImp ()
const {
return static_cast<const GeometryImp<mydim,cdim,GridImp>&
>(*this);}
350 template<
int cdim,
class Gr
idImp,
template<
int,
int,
class>
class GeometryImp>
361 typedef typename GridImp::ctype
ctype;
373 FieldVector<ctype, cdim>
global (
const FieldVector<ctype, mydim>&
local)
const 375 return asImp().corner(0);
379 FieldVector<ctype, mydim>
local (
const FieldVector<ctype, cdim>& )
const 381 return FieldVector<ctype, mydim>();
393 return asImp().corner(0);
398 GeometryImp<mydim,cdim,GridImp>& asImp () {
return static_cast<GeometryImp<mydim,cdim,GridImp>&
>(*this);}
399 const GeometryImp<mydim,cdim,GridImp>& asImp ()
const {
return static_cast<const GeometryImp<mydim,cdim,GridImp>&
>(*this);}
404 #endif // DUNE_GRID_GEOMETRY_HH FieldVector< ctype, cdim > GlobalCoordinate
Definition: common/geometry.hh:303
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:370
FieldVector< ctype, mydim > LocalCoordinate
Definition: common/geometry.hh:363
bool affine() const
Return true if the geometry mapping is affine and false otherwise.
Definition: common/geometry.hh:126
GridImp::ctype ctype
Definition: common/geometry.hh:300
Implementation::JacobianTransposed JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:118
FieldVector< ctype, cdim > center() const
return center of the geometry
Definition: common/geometry.hh:391
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
GeometryType type() const
Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.
Definition: common/geometry.hh:123
GlobalCoordinate corner(int i) const
Obtain a corner of the geometry.
Definition: common/geometry.hh:148
FieldVector< ctype, mydim > local(const FieldVector< ctype, cdim > &) const
return empty vector
Definition: common/geometry.hh:379
FieldVector< ctype, cdim > GlobalCoordinate
Definition: common/geometry.hh:364
FieldMatrix< ctype, mydim, cdim > JacobianTransposed
type of jacobian transposed
Definition: common/geometry.hh:309
Definition: common/geometry.hh:24
GlobalCoordinate center() const
return center of geometry
Definition: common/geometry.hh:215
LocalCoordinate local(const GlobalCoordinate &global) const
Evaluate the inverse map .
Definition: common/geometry.hh:166
GeometryImp< mydim, cdim, GridImp > Implementation
Definition: common/geometry.hh:78
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:367
const Implementation & impl() const
return reference to the implementation
Definition: common/geometry.hh:81
Default implementation for class Geometry.
Definition: common/geometry.hh:293
int corners() const
Return the number of corners of the reference element.
Definition: common/geometry.hh:134
Implementation realGeometry
Definition: common/geometry.hh:280
FieldMatrix< ctype, cdim, mydim > JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:306
Definition: common/geometry.hh:87
Include standard header files.
Definition: agrid.hh:59
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: common/geometry.hh:96
GridImp::ctype ctype
Definition: common/geometry.hh:361
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the map .
Definition: common/geometry.hh:157
Wrapper class for geometries.
Definition: common/geometry.hh:65
ctype integrationElement(const LocalCoordinate &local) const
Return the factor appearing in the integral transformation formula.
Definition: common/geometry.hh:194
FieldVector< ctype, cdim > global(const FieldVector< ctype, mydim > &local) const
return the only coordinate
Definition: common/geometry.hh:373
FieldVector< ctype, mydim > LocalCoordinate
Definition: common/geometry.hh:302
ctype volume() const
return volume of the geometry
Definition: common/geometry.hh:385
Implementation::JacobianInverseTransposed JacobianInverseTransposed
type of jacobian inverse transposed
Definition: common/geometry.hh:107
ctype volume() const
return volume of geometry
Definition: common/geometry.hh:200
GlobalCoordinate center() const
return center of the geometry
Definition: common/geometry.hh:331
ctype volume() const
return volume of the geometry
Definition: common/geometry.hh:312
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: common/geometry.hh:90
Geometry(const Implementation &impl)
copy constructor from implementation
Definition: common/geometry.hh:268
JacobianTransposed jacobianTransposed(const LocalCoordinate &local) const
Return the transposed of the Jacobian.
Definition: common/geometry.hh:231
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: common/geometry.hh:93
JacobianInverseTransposed jacobianInverseTransposed(const LocalCoordinate &local) const
Return inverse of transposed of Jacobian.
Definition: common/geometry.hh:257
Definition: common/geometry.hh:85