Dune Core Modules (2.10.0)
Wrapper class for geometries. More...
#include <dune/grid/common/geometry.hh>
Public Types | |
typedef GeometryImp< mydim, cdim, GridImp > | Implementation |
type of underlying implementation More... | |
typedef GridImp::ctype | ctype |
define type used for coordinates in grid module | |
typedef FieldVector< ctype, mydim > | LocalCoordinate |
type of local coordinates | |
typedef FieldVector< ctype, cdim > | GlobalCoordinate |
type of the global coordinates | |
typedef decltype(std::declval< Implementation >().volume()) | Volume |
Number type used for the geometry volume. | |
typedef Implementation::JacobianInverseTransposed | JacobianInverseTransposed |
type of jacobian inverse transposed More... | |
typedef Implementation::JacobianTransposed | JacobianTransposed |
type of jacobian transposed More... | |
using | JacobianInverse = Std::detected_or_t< JacobianInverseDefault, JacobianInverseOfImplementation, Implementation > |
type of jacobian inverse More... | |
using | Jacobian = Std::detected_or_t< JacobianDefault, JacobianOfImplementation, Implementation > |
type of jacobian More... | |
Public Member Functions | |
Implementation & | impl () |
access to the underlying implementation More... | |
const Implementation & | impl () 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. | |
bool | affine () const |
Return true if the geometry mapping is affine and false otherwise. | |
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 (const GlobalCoordinate &global) const | |
Evaluate the inverse map \( g^{-1}\). More... | |
Volume | integrationElement (const LocalCoordinate &local) const |
Return the factor appearing in the integral transformation formula. More... | |
Volume | volume () const |
return volume of geometry | |
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... | |
Jacobian | jacobian (const LocalCoordinate &local) const |
Return the Jacobian. More... | |
JacobianInverse | jacobianInverse (const LocalCoordinate &local) const |
Return inverse of Jacobian. More... | |
Static Public Attributes | |
static constexpr int | mydimension = mydim |
geometry dimension | |
static constexpr int | coorddimension = cdim |
dimension of embedding coordinate system | |
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 &) -> decltype(referenceElement< typename GridImp::ctype, mydim >(geo.type())) |
Second-level dispatch to select the correct reference element for a grid geometry. More... | |
Interface for grid implementers | |
Implementation | realGeometry |
Geometry (const Implementation &impl) | |
copy constructor from implementation | |
Detailed Description
class Dune::Geometry< mydim, cdim, GridImp, GeometryImp >
Wrapper class for geometries.
- Template Parameters
-
mydim Dimension of the domain cdim Dimension of the range GridImp Type that is a model of Dune::Grid GeometryImp Class 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 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.
@addtogroup Grid Grid The Dune Grid module defines a general interface to a parallel, in general
nonconforming, locally refined and hierarchical finite element mesh. The interface is independent of dimension and element type.
Terminology
@subsection subs1 Entity An entity is a geometric object that is part of a grid. It is
generalized polytope that has the same dimensionality as the grid or a lower dimension.
@subsection subs20 Dimension A grid has a fixed dimension \f$d\f$ which is the number of coordinates
required to specify any point in the grid. The dimension is a template parameter of a grid.
@subsection subs21 Codimension of an entity Each entity has a codimension \f$c\f$ where \f$0 \leq c \leq d\f$ (the dimension of the grid). An entity with codimension \f$ c\f$ in a grid of dimension \f$ d\f$ is a \f$d-c\f$-dimensional object. @subsection subs5 Subentity Entities are hierarchically constructed in the sense that entities of codimension 0 are made up of entities of codimension 1 which are themselves made up of entities of codimension 2 etc. until entities of codimension \f$d-1\f$ which consist of entities of codimension \f$ d\f$. @subsection subs3 Element An element is an entity of codimension 0. @subsection subs4 Vertex A vertex is an entity of codimension \f$ d\f$ (the same as the grid's dimension). @subsection subs22 World dimension Each grid has a world dimension \f$ w\f$ with \f$ w\geq d\f$. This is the number
of coordinates of the positions of the grid's vertices.
@subsection subs33 Hierarchical grid The %Dune grid interface describes not only a single grid but a sequence of grids with different resolution. This is achieved by beginning with an intentionally coarse grid, the so-called macro grid. Then each
element may be individually subdivided to yield new (smaller) elements. This construction is recursive such that each macro element and all the elements that resulted from subdividing it form a tree structure.
Grid refinement
The grid can only be modified in special phases, the so-called refinement phase. In between refinement phases the entities of the grid can not be modified in any way. During refinement currently only the hierachic subdivision can be modified.
@subsection subs3333 Grid level All elements of the macro grid form level 0 of the grid structure. All elements that are obtained from an \f$ l\f$-fold subdivision of a macro element form level \f$ l\f$ of the grid structure. @subsection subs333 Leaf grid All elements of a grid that are not subdivided any further make up the leaf grid. The leaf grid is the mesh with the finest resolution. @subsection subs6 Assignable A type is said to be assignable if it has a (public) copy constructor and assignment operator. Note that this definition requires always both methods. @subsection subs7 Default-constructible A type is said to be default-constructible if it has a constructor without arguments. @subsection subs8 Copy-constructible from type X A type is said to be copy constructible from some other type X if it has a copy constructor that takes a reference to an object of type X. @subsection subs9 Equality-comparable A type is said to be equality-comparable if it has an operator==. @subsection subs10 LessThan-comparable A type is lessthan-comparable if it has an operator<. @subsection subs11 Dereferenceable A type is dereferenceable if it has an operator* that delivers a reference to a value type. @subsection subs11 Iterator An iterator is a type that can be dereferenced to yield an object of its value type, i.e. it behaves like a pointer, and that can be incremented to point to the next element in a linear sequence. In that respect it is comparable to ForwardIterator in the Standard Template Library. @subsection subs12 Mutable iterator An iterator is called mutable if the value it refers to can be changed, i.e. it is assignable. @subsection subs13 Immutable iterator An iterator is called immutable if the value referenced by the iterator can not be changed, i. e. the value is not assignable and only methods marked const on the value can be called. @subsection subs14 Model A type M is called a model of another type X if it implements all the methods of X with the intended semantics. Typically X is a type that describes an interface. @section Grid3 Types common to all grid implementations - Dune::ReferenceElement describes the topology and geometry of standard entities. Any given entity of the grid can be completely specified by a reference element and a map from this reference element to world coordinate space. - Dune::GeometryType defines names for the reference elements. - Dune::Communication defines an interface to global communication operations in a portable and transparent way. In particular also for sequential grids. @section Grid2 Types making up a grid implementation
Each implementation of the Dune grid interface consist of a number of related types which together form a model of the grid interface. These types are the following:
- Grid which is a model of Dune::Grid where the template parameters are at least the dimension and the world dimension. It is a container of entities that allows to access these entities and that knows the number of entities.
- Entity which is a model of Dune::Entity. This class is parametrized by dimension and codimension. The entity encapsulates the topological part of an entity, i.e. its hierarchical construction from subentities and the relation to other entities. Entities cannot be created, copied or modified by the user. They can only be read-accessed through immutable iterators.
- Geometry which is a model of Dune::Geometry. This class encapsulates the geometric part of an entity by mapping local coordinates in a reference element to world coordinates.
- LevelIterator which is a model of Dune::LevelIterator is an immutable iterator that provides access to all entities of a given codimension and level of the grid.
- LeafIterator which is a model of Dune::LeafIterator is an immutable iterator that provides access to all entities of a given codimension of the leaf grid.
- HierarchicIterator which is a model of Dune::HierarchicIterator is an immutable iterator that provides access to all entities of codimension 0 that resulted from subdivision of a given entity of codimension 0.
- Intersection which is a model of Dune::Intersection provides access an intersection of codimension 1 of two entity of codimension 0 or one entity and the boundary. In a conforming mesh this is a face of an element. For two entities with a common intersection the Intersection also provides information about the geometric location of the intersection. Furthermore it also provides information about intersections of an entity with the internal or external boundaries.
- IntersectionIterator which is a model of Dune::IntersectionIterator provides access to all intersections of a given entity of codimension 0.
- LevelIndexSet and LeafIndexSet which are both models of Dune::IndexSet are used to attach any kind of user-defined data to (subsets of) entities of the grid. This data is supposed to be stored in one-dimensional arrays for reasons of efficiency.
- LocalIdSet and GlobalIdSet which are both models of Dune::IdSet are used to save user data during a grid refinement phase and during dynamic load balancing in the parallel case.
@section Grid22 Overview of basic capabilities of the types <TABLE> <TR> <TD>Class</TD> <TD>Assignable</TD> <TD>DefaultConstructible</TD> <TD>EqualityComparable</TD> <TD>LessThanComparable</TD> </TR> <TR> <TD>Grid</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> </TR> <TR> <TD>Entity</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> </TR> <TR> <TD>GeometryType</TD> <TD>yes</TD> <TD>yes</TD> <TD>yes</TD> <TD>yes</TD> </TR> <TR> <TD>Geometry</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> </TR> <TR> <TD>LevelIterator</TD> <TD>yes</TD> <TD>no</TD> <TD>yes</TD> <TD>no</TD> </TR> <TR> <TD>LeafIterator</TD> <TD>yes</TD> <TD>no</TD> <TD>yes</TD> <TD>no</TD> </TR> <TR> <TD>HierarchicIterator</TD> <TD>yes</TD> <TD>no</TD> <TD>yes</TD> <TD>no</TD> </TR> <TR> <TD>Intersection</TD> <TD>yes</TD> <TD>no</TD> <TD>yes</TD> <TD>no</TD> </TR> <TR> <TD>IntersectionIterator</TD> <TD>yes</TD> <TD>no</TD> <TD>yes</TD> <TD>no</TD> </TR> <TR> <TD>IndexSet</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> </TR> <TR> <TD>IdSet</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> <TD>no</TD> </TR> </TABLE>
Member Typedef Documentation
◆ Implementation
typedef GeometryImp< mydim, cdim, GridImp > Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::Implementation |
type of underlying implementation
- Warning
- Implementation details may change without prior notification.
◆ Jacobian
using Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::Jacobian = Std::detected_or_t<JacobianDefault, JacobianOfImplementation, Implementation> |
type of jacobian
The exact type is implementation-dependent. However, it is guaranteed to have the following properties:
- You can multiply it from the right to a suitable FieldMatrix
- It is copy constructible and copy assignable.
◆ JacobianInverse
using Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::JacobianInverse = Std::detected_or_t<JacobianInverseDefault, JacobianInverseOfImplementation, Implementation> |
type of jacobian inverse
The exact type is implementation-dependent. However, it is guaranteed to have the following properties:
- You can multiply it from the right to a suitable FieldMatrix
- It is copy constructible and copy assignable.
◆ JacobianInverseTransposed
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
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.
Member Function Documentation
◆ center()
|
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.
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
◆ corner()
|
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] i number of the corner (with respect to the reference element)
- Returns
- position of the i-th corner
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
Referenced by Dune::EdgeS0_5Basis< Geometry, RF >::EdgeS0_5Basis(), and Dune::EdgeS0_5Interpolation< Geometry, Traits_ >::EdgeS0_5Interpolation().
◆ corners()
|
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.
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
◆ global()
|
inline |
Evaluate the map \( g\).
- Parameters
-
[in] local Position in the reference element \(D\)
- Returns
- Position in \(W\)
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
Referenced by Dune::HierarchicSearch< Grid, IS >::findEntity(), and Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::LocalCoordinate().
◆ impl() [1/2]
|
inline |
access to the underlying implementation
- Warning
- Implementation details may change without prior notification.
Referenced by Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::affine(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::center(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::corner(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::corners(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::global(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::integrationElement(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::jacobianInverseTransposed(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::jacobianTransposed(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::LocalCoordinate(), Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::type(), and Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::volume().
◆ impl() [2/2]
|
inline |
access to the underlying implementation
- Warning
- Implementation details may change without prior notification.
◆ integrationElement()
|
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] local Position \(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}\).
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
◆ jacobian()
|
inline |
Return the Jacobian.
The Jacobian is defined in the documentation of integrationElement.
- Parameters
-
[in] local position \(x\in D\)
- Returns
- \(J_g(x)\)
- Note
- The exact return type is implementation defined.
◆ jacobianInverse()
|
inline |
Return inverse of Jacobian.
The Jacobian is defined in the documentation of integrationElement.
- Parameters
-
[in] local position \(x\in D\)
- Returns
- \(J_g^{-1}(x)\)
The use of this function is to compute the jacobians 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_{\hat{f}}(x) J_g^{-1}(x).\]
- Note
- In the non-quadratic case \(\textrm{cdim} \neq \textrm{mydim}\), the pseudoinverse of \(J_g(x)\) is returned. This means that its transposed is inverse for all tangential vectors in \(g(x)\) while mapping all normal vectors to zero.
- The exact return type is implementation defined.
◆ jacobianInverseTransposed()
|
inline |
Return inverse of transposed of Jacobian.
The Jacobian is defined in the documentation of integrationElement.
- Parameters
-
[in] local position \(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.
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
◆ jacobianTransposed()
|
inline |
Return the transposed of the Jacobian.
The Jacobian is defined in the documentation of integrationElement.
- Parameters
-
[in] local position \(x\in D\)
- Returns
- \(J_g^T(x)\)
- Note
- The exact return type is implementation defined.
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
◆ LocalCoordinate()
|
inline |
Evaluate the inverse map \( g^{-1}\).
- Parameters
-
[in] global Position in \(W\)
- Returns
- Position in \(D\) that maps to global
References Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::global(), and Dune::Geometry< mydim, cdim, GridImp, GeometryImp >::impl().
Friends And Related Function Documentation
◆ referenceElement()
|
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.
The documentation for this class was generated from the following file:
- dune/grid/common/geometry.hh