Dune Core Modules (2.9.0)

Dune::AffineGeometry< ct, mydim, cdim > Class Template Reference

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, mydimensionLocalCoordinate
 Type for local coordinate vector.
 
typedef FieldVector< ctype, coorddimensionGlobalCoordinate
 Type for coordinate vector in world space.
 
typedef ctype Volume
 Type used for volume.
 
typedef FieldMatrix< ctype, mydimension, coorddimensionJacobianTransposed
 Type for the transposed Jacobian matrix.
 
typedef FieldMatrix< ctype, coorddimension, mydimensionJacobianInverseTransposed
 Type for the transposed inverse Jacobian matrix.
 
typedef FieldMatrix< ctype, coorddimension, mydimensionJacobian
 Type for the Jacobian matrix.
 
typedef FieldMatrix< ctype, mydimension, coorddimensionJacobianInverse
 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 ([[maybe_unused]] const LocalCoordinate &local) const
 Obtain the integration element. More...
 
Volume volume () const
 Obtain the volume of the element.
 
const JacobianTransposedjacobianTransposed ([[maybe_unused]] const LocalCoordinate &local) const
 Obtain the transposed of the Jacobian. More...
 
const JacobianInverseTransposedjacobianInverseTransposed ([[maybe_unused]] const LocalCoordinate &local) const
 Obtain the transposed of the Jacobian's inverse. More...
 
Jacobian jacobian ([[maybe_unused]] const LocalCoordinate &local) const
 Obtain the Jacobian. More...
 
JacobianInverse jacobianInverse ([[maybe_unused]] 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

template<class ct, int mydim, int cdim>
class Dune::AffineGeometry< ct, mydim, cdim >

Implementation of the Geometry interface for affine geometries.

Template Parameters
ctType used for coordinates
mydimDimension of the geometry
cdimDimension of the world space

Member Function Documentation

◆ global()

template<class ct , int mydim, int cdim>
GlobalCoordinate Dune::AffineGeometry< ct, mydim, cdim >::global ( const LocalCoordinate local) const
inline

Evaluate the mapping.

Parameters
[in]locallocal coordinate to map
Returns
corresponding global coordinate

References Dune::DenseMatrix< MAT >::umtv().

Referenced by Dune::AffineGeometry< ct, mydim, cdim >::center(), Dune::AffineGeometry< ct, mydim, cdim >::corner(), and Dune::AffineGeometry< ct, mydim, cdim >::LocalCoordinate().

◆ integrationElement()

template<class ct , int mydim, int cdim>
ctype Dune::AffineGeometry< ct, mydim, cdim >::integrationElement ( [[maybe_unused] ] const LocalCoordinate local) const
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]locallocal coordinate to evaluate the integration element in
Returns
the integration element \(\mu(x)\).

◆ jacobian()

template<class ct , int mydim, int cdim>
Jacobian Dune::AffineGeometry< ct, mydim, cdim >::jacobian ( [[maybe_unused] ] const LocalCoordinate local) const
inline

Obtain the Jacobian.

Parameters
[in]locallocal coordinate to evaluate Jacobian in
Returns
a copy of the transposed of the Jacobian

References Dune::FieldMatrix< K, ROWS, COLS >::transposed().

◆ jacobianInverse()

template<class ct , int mydim, int cdim>
JacobianInverse Dune::AffineGeometry< ct, mydim, cdim >::jacobianInverse ( [[maybe_unused] ] const LocalCoordinate local) const
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()

template<class ct , int mydim, int cdim>
const JacobianInverseTransposed& Dune::AffineGeometry< ct, mydim, cdim >::jacobianInverseTransposed ( [[maybe_unused] ] const LocalCoordinate local) const
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()

template<class ct , int mydim, int cdim>
const JacobianTransposed& Dune::AffineGeometry< ct, mydim, cdim >::jacobianTransposed ( [[maybe_unused] ] const LocalCoordinate local) const
inline

Obtain the transposed of the Jacobian.

Parameters
[in]locallocal coordinate to evaluate Jacobian in
Returns
a reference to the transposed of the Jacobian

◆ LocalCoordinate()

template<class ct , int mydim, int cdim>
Dune::AffineGeometry< ct, mydim, cdim >::LocalCoordinate ( const GlobalCoordinate global) const
inline

Evaluate the inverse mapping.

Parameters
[in]globalglobal coordinate to map
Returns
corresponding local coordinate

The returned local coordinate y minimizes

(global( y ) - x).two_norm()
GlobalCoordinate global(const LocalCoordinate &local) const
Evaluate the mapping.
Definition: affinegeometry.hh:581

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:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (Mar 28, 23:30, 2024)