Dune Core Modules (2.9.1)
Implementation of the Geometry interface for affine geometries. More...
#include <dune/geometry/affinegeometry.hh>
Public Types | |
typedef ct | ctype |
Type used for coordinates. | |
typedef FieldVector< ctype, mydimension > | LocalCoordinate |
Type for local coordinate vector. | |
typedef FieldVector< ctype, coorddimension > | GlobalCoordinate |
Type for coordinate vector in world space. | |
typedef ctype | Volume |
Type used for volume. | |
typedef FieldMatrix< ctype, mydimension, coorddimension > | JacobianTransposed |
Type for the transposed Jacobian matrix. | |
typedef FieldMatrix< ctype, coorddimension, mydimension > | JacobianInverseTransposed |
Type for the transposed inverse Jacobian matrix. | |
typedef FieldMatrix< ctype, coorddimension, mydimension > | Jacobian |
Type for the Jacobian matrix. | |
typedef FieldMatrix< ctype, mydimension, coorddimension > | JacobianInverse |
Type for the inverse Jacobian matrix. | |
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. | |
AffineGeometry (Dune::GeometryType gt, const GlobalCoordinate &origin, const JacobianTransposed &jt) | |
Create affine geometry from GeometryType, one vertex, and the Jacobian matrix. | |
template<class CoordVector > | |
AffineGeometry (const ReferenceElement &refElement, const CoordVector &coordVector) | |
Create affine geometry from reference element and a vector of vertex coordinates. | |
template<class CoordVector > | |
AffineGeometry (Dune::GeometryType gt, const CoordVector &coordVector) | |
Create affine geometry from GeometryType and a vector of vertex coordinates. | |
bool | affine () const |
Always true: this is an affine geometry. | |
Dune::GeometryType | type () const |
Obtain the type of the reference element. | |
int | corners () const |
Obtain number of corners of the corresponding reference element. | |
GlobalCoordinate | corner (int i) const |
Obtain coordinates of the i-th corner. | |
GlobalCoordinate | center () const |
Obtain the centroid of the mapping's image. | |
GlobalCoordinate | global (const LocalCoordinate &local) const |
Evaluate the mapping. More... | |
LocalCoordinate (const GlobalCoordinate &global) const | |
Evaluate the inverse mapping. More... | |
ctype | integrationElement (const LocalCoordinate &local) const |
Obtain the integration element. More... | |
Volume | volume () const |
Obtain the volume of the element. | |
const JacobianTransposed & | jacobianTransposed (const LocalCoordinate &local) const |
Obtain the transposed of the Jacobian. More... | |
const JacobianInverseTransposed & | jacobianInverseTransposed (const LocalCoordinate &local) const |
Obtain the transposed of the Jacobian's inverse. More... | |
Jacobian | jacobian (const LocalCoordinate &local) const |
Obtain the Jacobian. More... | |
JacobianInverse | jacobianInverse (const LocalCoordinate &local) const |
Obtain the Jacobian's inverse. More... | |
Static Public Attributes | |
static const int | mydimension = mydim |
Dimension of the geometry. | |
static const int | coorddimension = cdim |
Dimension of the world space. | |
Detailed Description
class Dune::AffineGeometry< ct, mydim, cdim >
Implementation of the Geometry interface for affine geometries.
- Template Parameters
-
ct Type used for coordinates mydim Dimension of the geometry cdim Dimension of the world space
Member Function Documentation
◆ global()
|
inline |
Evaluate the mapping.
- Parameters
-
[in] local local coordinate to map
- Returns
- corresponding global coordinate
References Dune::AffineGeometry< ct, mydim, cdim >::global(), and Dune::DenseMatrix< MAT >::umtv().
Referenced by Dune::AffineGeometry< ct, mydim, cdim >::center(), Dune::AffineGeometry< ct, mydim, cdim >::corner(), Dune::AffineGeometry< ct, mydim, cdim >::global(), and Dune::AffineGeometry< ct, mydim, cdim >::LocalCoordinate().
◆ integrationElement()
|
inline |
Obtain the integration element.
If the Jacobian of the mapping is denoted by $J(x)$, the integration integration element \(\mu(x)\) is given by
\[ \mu(x) = \sqrt{|\det (J^T(x) J(x))|}.\]
- Parameters
-
[in] local local coordinate to evaluate the integration element in
- Returns
- the integration element \(\mu(x)\).
◆ jacobian()
|
inline |
Obtain the Jacobian.
- Parameters
-
[in] local local coordinate to evaluate Jacobian in
- Returns
- a copy of the transposed of the Jacobian
References Dune::FieldMatrix< K, ROWS, COLS >::transposed().
◆ jacobianInverse()
|
inline |
Obtain the Jacobian's inverse.
The Jacobian's inverse is defined as a pseudo-inverse. If we denote the Jacobian by \(J(x)\), the following condition holds:
\[J^{-1}(x) J(x) = I.\]
References Dune::FieldMatrix< K, ROWS, COLS >::transposed().
◆ jacobianInverseTransposed()
|
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 \(J(x)\), the following condition holds:
\[J^{-1}(x) J(x) = I.\]
◆ jacobianTransposed()
|
inline |
Obtain the transposed of the Jacobian.
- Parameters
-
[in] local local coordinate to evaluate Jacobian in
- Returns
- a reference to the transposed of the Jacobian
◆ LocalCoordinate()
|
inline |
Evaluate the inverse mapping.
- Parameters
-
[in] global global coordinate to map
- Returns
- corresponding local coordinate
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.
References Dune::AffineGeometry< ct, mydim, cdim >::global(), and Dune::DenseMatrix< MAT >::mtv().
The documentation for this class was generated from the following file:
- dune/geometry/affinegeometry.hh