Dune Core Modules (2.10.0)
generic geometry implementation based on corner coordinates More...
#include <dune/geometry/multilineargeometry.hh>
Public Types | |
typedef ct | ctype |
coordinate type | |
typedef FieldVector< ctype, mydimension > | LocalCoordinate |
type of local coordinates | |
typedef FieldVector< ctype, coorddimension > | GlobalCoordinate |
type of global coordinates | |
typedef ctype | Volume |
type of volume | |
typedef FieldMatrix< ctype, mydimension, coorddimension > | JacobianTransposed |
type of jacobian transposed | |
typedef FieldMatrix< ctype, coorddimension, mydimension > | Jacobian |
Type for the Jacobian matrix. | |
typedef FieldMatrix< ctype, mydimension, coorddimension > | JacobianInverse |
Type for the inverse Jacobian matrix. | |
typedef ReferenceElements::ReferenceElement | ReferenceElement |
type of reference element | |
Public Member Functions | |
template<class Corners > | |
MultiLinearGeometry (const ReferenceElement &refElement, const Corners &corners) | |
constructor More... | |
template<class Corners > | |
MultiLinearGeometry (Dune::GeometryType gt, const Corners &corners) | |
constructor More... | |
bool | affine () const |
is this mapping affine? | |
Dune::GeometryType | type () const |
obtain the name 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 &globalCoord) const | |
evaluate the inverse mapping More... | |
Volume | integrationElement (const LocalCoordinate &local) const |
obtain the integration element More... | |
Volume | volume () const |
obtain the volume of the mapping's image More... | |
JacobianTransposed | jacobianTransposed (const LocalCoordinate &local) const |
obtain the transposed of the Jacobian More... | |
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 |
geometry dimension | |
static const int | coorddimension = cdim |
coordinate dimension | |
Detailed Description
class Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >
generic geometry implementation based on corner coordinates
Based on the recursive definition of the reference elements, the MultiLinearGeometry provides a generic implementation of a geometry given the corner coordinates.
The geometric mapping is multilinear in the classical sense only in the case of cubes; for simplices it is linear. The name is still justified, because the mapping satisfies the important property of begin linear along edges.
- Template Parameters
-
ct coordinate type mydim geometry dimension cdim coordinate dimension Traits traits allowing to tweak some implementation details (optional)
The requirements on the traits are documented along with their default, MultiLinearGeometryTraits.
Constructor & Destructor Documentation
◆ MultiLinearGeometry() [1/2]
|
inline |
constructor
- Parameters
-
[in] refElement reference element for the geometry [in] corners corners to store internally
- Note
- The type of corners is actually a template argument. It is only required that the internal corner storage can be constructed from this object.
◆ MultiLinearGeometry() [2/2]
|
inline |
constructor
- Parameters
-
[in] gt geometry type [in] corners corners to store internally
- Note
- The type of corners is actually a template argument. It is only required that the internal corner storage can be constructed from this object.
Member Function Documentation
◆ global()
|
inline |
evaluate the mapping
- Parameters
-
[in] local local coordinate to map
- Returns
- corresponding global coordinate
References Dune::get().
Referenced by Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::center(), Dune::P1VTKFunction< GV, V >::evaluate(), and Dune::CachedMultiLinearGeometry< ct, mydim, cdim, Traits >::global().
◆ 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)\).
- Note
- For affine mappings, it is more efficient to call jacobianInverseTransposed before integrationElement, if both are required.
References Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianTransposed().
Referenced by Dune::CachedMultiLinearGeometry< ct, mydim, cdim, Traits >::integrationElement(), and Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::volume().
◆ 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::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianTransposed(), and 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::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianInverseTransposed().
◆ 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.\]
Referenced by Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianInverse(), and Dune::CachedMultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianInverseTransposed().
◆ 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
- Note
- The returned reference is reused on the next call to JacobianTransposed, destroying the previous value.
References Dune::get().
Referenced by Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::integrationElement(), Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobian(), Dune::CachedMultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianTransposed(), and Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::LocalCoordinate().
◆ LocalCoordinate()
|
inline |
evaluate the inverse mapping
- Parameters
-
[in] globalCoord global coordinate to map
- Returns
- corresponding local coordinate
- Note
- For given global coordinate y the returned local coordinate x that minimizes the following function over the local coordinate space spanned by the reference element. (global( x ) - y).two_norm()GlobalCoordinate global(const LocalCoordinate &local) constevaluate the mappingDefinition: multilineargeometry.hh:290
References Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::affine(), Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::jacobianTransposed(), Dune::Hybrid::max, and Dune::DenseVector< V >::two_norm2().
◆ volume()
|
inline |
obtain the volume of the mapping's image
- Note
- The current implementation just returns which is wrong for n-linear surface maps and other nonlinear maps.integrationElement( refElement().position( 0, 0 ) ) * refElement().volume()Volume integrationElement(const LocalCoordinate &local) constobtain the integration elementDefinition: multilineargeometry.hh:350
References Dune::MultiLinearGeometry< ct, mydim, cdim, Traits >::integrationElement().
Referenced by Dune::CachedMultiLinearGeometry< ct, mydim, cdim, Traits >::volume().
The documentation for this class was generated from the following file:
- dune/geometry/multilineargeometry.hh