dune-grid  2.6-git
Public Types | Public Member Functions | Protected Attributes | Related Functions | List of all members
Dune::Geometry< mydim, cdim, GridImp, GeometryImp > Class Template Reference

Wrapper class for geometries. More...

#include <dune/grid/common/geometry.hh>

Public Types

enum  { mydimension =mydim }
 export geometry dimension More...
 
enum  { coorddimension =cdim }
 export coordinate dimension More...
 
typedef GeometryImp< mydim, cdim, GridImp > Implementation
 type of underlying implementation More...
 
typedef GridImp::ctype ctype
 define type used for coordinates in grid module More...
 
typedef FieldVector< ctype, mydim > LocalCoordinate
 type of local coordinates More...
 
typedef FieldVector< ctype, cdim > GlobalCoordinate
 type of the global coordinates More...
 
typedef Implementation::JacobianInverseTransposed JacobianInverseTransposed
 type of jacobian inverse transposed More...
 
typedef Implementation::JacobianTransposed JacobianTransposed
 type of jacobian transposed More...
 

Public Member Functions

Implementationimpl ()
 access to the underlying implementation More...
 
const Implementationimpl () const
 access to the underlying implementation More...
 
GeometryType type () const
 Return the type of the reference element. The type can be used to access the Dune::ReferenceElement. More...
 
bool affine () const
 Return true if the geometry mapping is affine and false otherwise. More...
 
int corners () const
 Return the number of corners of the reference element. More...
 
GlobalCoordinate corner (int i) const
 Obtain a corner of the geometry. More...
 
GlobalCoordinate global (const LocalCoordinate &local) const
 Evaluate the map $ g$. More...
 
LocalCoordinate local (const GlobalCoordinate &global) const
 Evaluate the inverse map $ g^{-1}$. More...
 
ctype integrationElement (const LocalCoordinate &local) const
 Return the factor appearing in the integral transformation formula. More...
 
ctype volume () const
 return volume of geometry More...
 
GlobalCoordinate center () const
 return center of geometry More...
 
JacobianTransposed jacobianTransposed (const LocalCoordinate &local) const
 Return the transposed of the Jacobian. More...
 
JacobianInverseTransposed jacobianInverseTransposed (const LocalCoordinate &local) const
 Return inverse of transposed of Jacobian. More...
 
Interface for grid implementers
 Geometry (const Implementation &impl)
 copy constructor from implementation More...
 

Protected Attributes

Implementation realGeometry
 

Related Functions

(Note that these are not member functions.)

template<int mydim, int cdim, class GridImp , template< int, int, class > class GeometryImp, typename Impl >
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. More...
 

Detailed Description

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
class Dune::Geometry< mydim, cdim, GridImp, GeometryImp >

Wrapper class for geometries.

Template Parameters
mydimDimension of the domain
cdimDimension of the range
GridImpType that is a model of Dune::Grid
GeometryImpClass template that is a model of Dune::Geometry

Maps

A Geometry defines a map

\[ g : D \to W\]

where $D\subseteq\mathbf{R}^\textrm{mydim}$ and $W\subseteq\mathbf{R}^\textrm{cdim}$. The domain $D$ is one of a set of predefined convex polytopes, the so-called reference elements (

See also
Dune::ReferenceElement). The dimensionality of $D$ is mydim. In general $\textrm{mydim}\leq\textrm{cdim}$, i.e. the convex polytope may be mapped to a manifold. Moreover, we require that $ g\in \left( C^1(D) \right)^\textrm{cdim}$ and one-to-one.

Engine Concept

The Geometry class template wraps an object of type GeometryImp and forwards all member function calls to corresponding members of this class. In that sense Geometry defines the interface and GeometryImp supplies the implementation.

Member Typedef Documentation

◆ ctype

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef GridImp::ctype Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::ctype

define type used for coordinates in grid module

◆ GlobalCoordinate

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef FieldVector< ctype, cdim > Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::GlobalCoordinate

type of the global coordinates

◆ Implementation

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef GeometryImp< mydim, cdim, GridImp > Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::Implementation

type of underlying implementation

Warning
Implementation details may change without prior notification.

◆ JacobianInverseTransposed

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef Implementation::JacobianInverseTransposed Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::JacobianInverseTransposed

type of jacobian inverse transposed

The exact type is implementation-dependent. However, it is guaranteed to have the following properties:

  • It satisfies the ConstMatrix interface.
  • It is copy constructible and copy assignable.

◆ JacobianTransposed

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef Implementation::JacobianTransposed Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::JacobianTransposed

type of jacobian transposed

The exact type is implementation-dependent. However, it is guaranteed to have the following properties:

  • It satisfies the ConstMatrix interface.
  • It is copy constructible and copy assignable.

◆ LocalCoordinate

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
typedef FieldVector<ctype, mydim> Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::LocalCoordinate

type of local coordinates

Member Enumeration Documentation

◆ anonymous enum

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
anonymous enum

export geometry dimension

Enumerator
mydimension 

geometry dimension

◆ anonymous enum

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
anonymous enum

export coordinate dimension

Enumerator
coorddimension 

dimension of embedding coordinate system

Constructor & Destructor Documentation

◆ Geometry()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::Geometry ( const Implementation impl)
inlineexplicit

copy constructor from implementation

Member Function Documentation

◆ affine()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
bool Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::affine ( ) const
inline

Return true if the geometry mapping is affine and false otherwise.

◆ center()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
GlobalCoordinate Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::center ( ) const
inline

return center of geometry

Note that this method is still subject to a change of name and semantics. At the moment, the center is not required to be the centroid of the geometry, or even the centroid of its corners. This makes the current default implementation acceptable, which maps the centroid of the reference element to the geometry. We may change the name (and semantic) of the method to centroid() if we find reasonably efficient ways to implement it properly.

◆ corner()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
GlobalCoordinate Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::corner ( int  i) const
inline

Obtain a corner of the geometry.

This method is for convenient access to the corners of the geometry. The same result could be achieved by calling

Parameters
[in]inumber of the corner (with respect to the reference element)
Returns
position of the i-th corner

◆ corners()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
int Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::corners ( ) const
inline

Return the number of corners of the reference element.

Since a geometry is a convex polytope the number of corners is a well-defined concept. The method is redundant because this information is also available via the reference element. It is here for efficiency and ease of use.

◆ global()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
GlobalCoordinate Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::global ( const LocalCoordinate local) const
inline

Evaluate the map $ g$.

Parameters
[in]localPosition in the reference element $D$
Returns
Position in $W$

◆ impl() [1/2]

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
Implementation& Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl ( )
inline

access to the underlying implementation

Warning
Implementation details may change without prior notification.

◆ impl() [2/2]

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
const Implementation& Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl ( ) const
inline

access to the underlying implementation

Warning
Implementation details may change without prior notification.

◆ integrationElement()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
ctype Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::integrationElement ( const LocalCoordinate local) const
inline

Return the factor appearing in the integral transformation formula.

Let $ g : D \to W$ denote the transformation described by the Geometry. Then the jacobian of the transformation is defined as the $\textrm{cdim}\times\textrm{mydim}$ matrix

\[ J_g(x) = \left( \begin{array}{ccc} \frac{\partial g_0}{\partial x_0} & \cdots & \frac{\partial g_0}{\partial x_{n-1}} \\ \vdots & \ddots & \vdots \\ \frac{\partial g_{m-1}}{\partial x_0} & \cdots & \frac{\partial g_{m-1}}{\partial x_{n-1}} \end{array} \right).\]

Here we abbreviated $m=\textrm{cdim}$ and $n=\textrm{mydim}$ for ease of readability.

The integration element $\mu(x)$ for any $x\in D$ is then defined as

\[ \mu(x) = \sqrt{|\det J_g^T(x)J_g(x)|}.\]

Parameters
[in]localPosition $x\in D$
Returns
integration element $\mu(x)$
Note
Each implementation computes the integration element with optimal efficiency. For example in an equidistant structured mesh it may be as simple as $h^\textrm{mydim}$.

◆ jacobianInverseTransposed()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
JacobianInverseTransposed Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::jacobianInverseTransposed ( const LocalCoordinate local) const
inline

Return inverse of transposed of Jacobian.

The Jacobian is defined in the documentation of integrationElement.

Parameters
[in]localposition $x\in D$
Returns
$J_g^{-T}(x)$

The use of this function is to compute the gradient of some function $f : W \to \textbf{R}$ at some position $y=g(x)$, where $x\in D$ and $g$ the transformation of the Geometry. When we set $\hat{f}(x) = f(g(x))$ and apply the chain rule we obtain

\[\nabla f(g(x)) = J_g^{-T}(x) \nabla \hat{f}(x).\]

Note
In the non-quadratic case $\textrm{cdim} \neq \textrm{mydim}$, the pseudoinverse of $J_g^T(x)$ is returned. This means that it is inverse for all tangential vectors in $g(x)$ while mapping all normal vectors to zero.
The exact return type is implementation defined.

◆ jacobianTransposed()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
JacobianTransposed Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::jacobianTransposed ( const LocalCoordinate local) const
inline

Return the transposed of the Jacobian.

The Jacobian is defined in the documentation of integrationElement.

Parameters
[in]localposition $x\in D$
Returns
$J_g^T(x)$
Note
The exact return type is implementation defined.

◆ local()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
LocalCoordinate Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::local ( const GlobalCoordinate global) const
inline

Evaluate the inverse map $ g^{-1}$.

Parameters
[in]globalPosition in $W$
Returns
Position in $D$ that maps to global

◆ type()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
GeometryType Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::type ( ) const
inline

Return the type of the reference element. The type can be used to access the Dune::ReferenceElement.

◆ volume()

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
ctype Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::volume ( ) const
inline

return volume of geometry

Friends And Related Function Documentation

◆ referenceElement()

template<int mydim, int cdim, class GridImp , template< int, int, class > class GeometryImp, typename Impl >
auto referenceElement ( const Geometry< mydim, cdim, GridImp, GeometryImp > &  geo,
const Impl &  impl 
) -> decltype(referenceElement<typename GridImp::ctype,mydim>(geo.type()))
related

Second-level dispatch to select the correct reference element for a grid geometry.

This function is the default implementation of the second-level reference element dispatch performed by Geometry.

Note
This function is only important for grid implementors with geometries that do not have a standard reference element.

When referenceElement() is called with a Geometry, it will forward the call to referenceElement(const Geometry&,const GeometryImplementation&). This default implementation will do the right thing as long as your geometry is based on a standard Dune ReferenceElement. If it is not and you want to supply your own reference element implementation, provide an override of this function for your specific geometry implementation.

Member Data Documentation

◆ realGeometry

template<int mydim, int cdim, class GridImp, template< int, int, class > class GeometryImp>
Implementation Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::realGeometry
protected

The documentation for this class was generated from the following file: