|
dune-geometry
2.8.0
|
Implementation of the Geometry interface for affine geometries. More...
#include <dune/geometry/affinegeometry.hh>
Public Types | |
| typedef ct | ctype |
| Type used for coordinates. More... | |
| typedef FieldVector< ctype, mydimension > | LocalCoordinate |
| Type for local coordinate vector. More... | |
| typedef FieldVector< ctype, coorddimension > | GlobalCoordinate |
| Type for coordinate vector in world space. More... | |
| typedef ctype | Volume |
| Type used for volume. More... | |
| typedef FieldMatrix< ctype, mydimension, coorddimension > | JacobianTransposed |
| Type for the transposed Jacobian matrix. More... | |
| typedef FieldMatrix< ctype, coorddimension, mydimension > | JacobianInverseTransposed |
| Type for the transposed inverse Jacobian matrix. More... | |
Public Member Functions | |
| AffineGeometry (const ReferenceElement &refElement, const GlobalCoordinate &origin, const JacobianTransposed &jt) | |
| Create affine geometry from reference element, one vertex, and the Jacobian matrix. More... | |
| AffineGeometry (Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt) | |
| Create affine geometry from GeometryType, one vertex, and the Jacobian matrix. More... | |
| template<class CoordVector > | |
| AffineGeometry (const ReferenceElement &refElement, const CoordVector &coordVector) | |
| Create affine geometry from reference element and a vector of vertex coordinates. More... | |
| template<class CoordVector > | |
| AffineGeometry (Dune::GeometryType gt, const CoordVector &coordVector) | |
| Create affine geometry from GeometryType and a vector of vertex coordinates. More... | |
| bool | affine () const |
| Always true: this is an affine geometry. More... | |
| Dune::GeometryType | type () const |
| Obtain the type of the reference element. More... | |
| int | corners () const |
| Obtain number of corners of the corresponding reference element. More... | |
| GlobalCoordinate | corner (int i) const |
| Obtain coordinates of the i-th corner. More... | |
| GlobalCoordinate | center () const |
| Obtain the centroid of the mapping's image. More... | |
| GlobalCoordinate | global (const LocalCoordinate &local) const |
| Evaluate the mapping. More... | |
| LocalCoordinate (const GlobalCoordinate &global) const | |
| Evaluate the inverse mapping. More... | |
| ctype | integrationElement ([[maybe_unused]] const LocalCoordinate &local) const |
| Obtain the integration element. More... | |
| Volume | volume () const |
| Obtain the volume of the element. More... | |
| const JacobianTransposed & | jacobianTransposed ([[maybe_unused]] const LocalCoordinate &local) const |
| Obtain the transposed of the Jacobian. More... | |
| const JacobianInverseTransposed & | jacobianInverseTransposed ([[maybe_unused]] const LocalCoordinate &local) const |
| Obtain the transposed of the Jacobian's inverse. More... | |
Static Public Attributes | |
| static const int | mydimension = mydim |
| Dimension of the geometry. More... | |
| static const int | coorddimension = cdim |
| Dimension of the world space. More... | |
Implementation of the Geometry interface for affine geometries.
| ct | Type used for coordinates |
| mydim | Dimension of the geometry |
| cdim | Dimension of the world space |
| typedef ct Dune::AffineGeometry< ct, mydim, cdim >::ctype |
Type used for coordinates.
| typedef FieldVector< ctype, coorddimension > Dune::AffineGeometry< ct, mydim, cdim >::GlobalCoordinate |
Type for coordinate vector in world space.
| typedef FieldMatrix< ctype, coorddimension, mydimension > Dune::AffineGeometry< ct, mydim, cdim >::JacobianInverseTransposed |
Type for the transposed inverse Jacobian matrix.
| typedef FieldMatrix< ctype, mydimension, coorddimension > Dune::AffineGeometry< ct, mydim, cdim >::JacobianTransposed |
Type for the transposed Jacobian matrix.
| typedef FieldVector< ctype, mydimension > Dune::AffineGeometry< ct, mydim, cdim >::LocalCoordinate |
Type for local coordinate vector.
| typedef ctype Dune::AffineGeometry< ct, mydim, cdim >::Volume |
Type used for volume.
|
inline |
Create affine geometry from reference element, one vertex, and the Jacobian matrix.
|
inline |
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix.
|
inline |
Create affine geometry from reference element and a vector of vertex coordinates.
|
inline |
Create affine geometry from GeometryType and a vector of vertex coordinates.
|
inline |
Always true: this is an affine geometry.
|
inline |
Obtain the centroid of the mapping's image.
|
inline |
Obtain coordinates of the i-th corner.
|
inline |
Obtain number of corners of the corresponding reference element.
|
inline |
Evaluate the mapping.
| [in] | local | local coordinate to map |
|
inline |
Obtain the integration element.
If the Jacobian of the mapping is denoted by $J(x)$, the integration integration element
is given by
| [in] | local | local coordinate to evaluate the integration element in |
|
inline |
Obtain the transposed of the Jacobian's inverse.
The Jacobian's inverse is defined as a pseudo-inverse. If we denote the Jacobian by
, the following condition holds:
|
inline |
Obtain the transposed of the Jacobian.
| [in] | local | local coordinate to evaluate Jacobian in |
|
inline |
Evaluate the inverse mapping.
| [in] | global | global coordinate to map |
The returned local coordinate y minimizes
on the entire affine hull of the reference element. This degenerates to the inverse map if the argument y is in the range of the map.
|
inline |
Obtain the type of the reference element.
|
inline |
Obtain the volume of the element.
|
static |
Dimension of the world space.
|
static |
Dimension of the geometry.