Dune Core Modules (unstable)

Dune::MappedGeometry< Map, Geo > Class Template Reference

Geometry parametrized by a LocalFunction and a LocalGeometry. More...

#include <dune/geometry/mappedgeometry.hh>

Public Types

using LocalCoordinate = typename Geo::LocalCoordinate
 type of local coordinates
 
using GlobalCoordinate = std::remove_reference_t< decltype(std::declval< Map >()(std::declval< typename Geo::GlobalCoordinate >()))>
 type of global coordinates
 
using ctype = typename Geo::ctype
 coordinate type
 
using Volume = std::remove_reference_t< decltype(Dune::power(std::declval< ctype >(), mydimension))>
 type of volume
 
using Jacobian = FieldMatrix< ctype, coorddimension, mydimension >
 type of jacobian
 
using JacobianTransposed = FieldMatrix< ctype, mydimension, coorddimension >
 type of jacobian transposed
 
using JacobianInverse = FieldMatrix< ctype, mydimension, coorddimension >
 type of jacobian inverse
 
using JacobianInverseTransposed = FieldMatrix< ctype, coorddimension, mydimension >
 type of jacobian inverse transposed
 

Public Member Functions

template<class Geo_ , class Map_ , std::enable_if_t< Dune::IsInteroperable< Map, Map_ >::value, int > = 0, std::enable_if_t< Dune::IsInteroperable< Geo, Geo_ >::value, int > = 0>
 MappedGeometry (Map_ &&mapping, Geo_ &&geometry, bool affine=false)
 Constructor from mapping to parametrize the geometry. More...
 
bool affine () const
 Is this mapping affine? Not in general, since we don't know anything about the mapping. The returned value can be given in the constructor of the class.
 
GeometryType type () const
 Obtain the geometry type from 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
 Map the center of the wrapped geometry.
 
GlobalCoordinate global (const LocalCoordinate &local) const
 Evaluate the coordinate mapping. More...
 
 LocalCoordinate (const GlobalCoordinate &y, Impl::GaussNewtonOptions< ctype > opts={}) const
 Evaluate the inverse coordinate mapping. More...
 
ctype integrationElement (const LocalCoordinate &local) const
 Obtain the integration element. More...
 
Volume volume (Impl::ConvergenceOptions< ctype > opts={}) const
 Obtain the volume of the mapping's image. More...
 
Volume volume (const QuadratureRule< ctype, mydimension > &quadRule) const
 Obtain the volume of the mapping's image by given quadrature rules.
 
Jacobian jacobian (const LocalCoordinate &local) const
 Obtain the Jacobian. More...
 
JacobianTransposed jacobianTransposed (const LocalCoordinate &local) const
 Obtain the transposed of the Jacobian. More...
 
JacobianInverse jacobianInverse (const LocalCoordinate &local) const
 Obtain the Jacobian's inverse. More...
 
JacobianInverseTransposed jacobianInverseTransposed (const LocalCoordinate &local) const
 Obtain the transposed of the Jacobian's inverse. More...
 

Static Public Attributes

static constexpr int mydimension = LocalCoordinate::size()
 geometry dimension
 
static constexpr int coorddimension = GlobalCoordinate::size()
 coordinate dimension
 

Detailed Description

template<class Map, class Geo>
class Dune::MappedGeometry< Map, Geo >

Geometry parametrized by a LocalFunction and a LocalGeometry.

This class represents a geometry that is parametrized by the chained mapping of a geometry g of type Geo and the given (differentiable) function f of type Map as f o g.

The geometry g: LG -> GG maps points from the local domain LG to its global domain GG. It must fulfill the dune-grid geometry-concept. Examples are the geometry of a grid element, the geometry of a `ReferenceElement` or any other geometry implementation in dune-geometry. The local domain LG represents the local coordinates of the MappedGeometry.

The function f: LF -> GF is a differentiable function in the sense of the dune-functions concept. It maps coordinates from the global domain of the geometry, LF=GG into the range type GF representing the global coordinates of the MappedGeometry.

Requirements: Let local be of type LG, g of type Geo and f of type Map.

  • g must be a model of Concept::Geometry
  • f must be a model of Concept::DifferentiableFunction<GF(LF)>

The following expressions must be valid:

  • GG gg = g.global(local): the "evaluation" of the geometry in its local domain.
  • GF gf = f(gg): the evaluation of the mapping f in its local domain that is the geometries' range.
  • auto df = derivative(f), df(gg): The derivative of f w.r.t. its global coordinate GF evaluated at its local coordinated LF=GG.
Template Parameters
MapDifferentiable mapping f: LF -> GF with LF = GG the range of the geometry g.
GeoA geometry type g: LG -> GG fulfilling the concept Concept::Geometry.

Constructor & Destructor Documentation

◆ MappedGeometry()

template<class Map , class Geo >
template<class Geo_ , class Map_ , std::enable_if_t< Dune::IsInteroperable< Map, Map_ >::value, int > = 0, std::enable_if_t< Dune::IsInteroperable< Geo, Geo_ >::value, int > = 0>
Dune::MappedGeometry< Map, Geo >::MappedGeometry ( Map_ &&  mapping,
Geo_ &&  geometry,
bool  affine = false 
)
inline

Constructor from mapping to parametrize the geometry.

Parameters
[in]mappingA differentiable function for the parametrization of the geometry
[in]geometryThe geometry that is mapped
[in]affineFlag indicating whether the MappedGeometry represents an affine mapping

Member Function Documentation

◆ global()

template<class Map , class Geo >
GlobalCoordinate Dune::MappedGeometry< Map, Geo >::global ( const LocalCoordinate local) const
inline

Evaluate the coordinate mapping.

Chained evaluation of the geometry and mapping from local coordinates in the reference element associated to the wrapped geometry.

Parameters
[in]locallocal coordinates
Returns
corresponding global coordinate

◆ integrationElement()

template<class Map , class Geo >
ctype Dune::MappedGeometry< Map, Geo >::integrationElement ( const LocalCoordinate local) const
inline

Obtain the integration element.

If the Jacobian of the geometry is denoted by \(J(x)\), the 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)\).

References Dune::MappedGeometry< Map, Geo >::jacobianTransposed().

Referenced by Dune::MappedGeometry< Map, Geo >::volume().

◆ jacobian()

template<class Map , class Geo >
Jacobian Dune::MappedGeometry< Map, Geo >::jacobian ( const LocalCoordinate local) const
inline

Obtain the Jacobian.

Parameters
[in]locallocal coordinate to evaluate Jacobian in
Returns
the matrix corresponding to the Jacobian

Referenced by Dune::MappedGeometry< Map, Geo >::jacobianInverse(), and Dune::MappedGeometry< Map, Geo >::jacobianTransposed().

◆ jacobianInverse()

template<class Map , class Geo >
JacobianInverse Dune::MappedGeometry< Map, Geo >::jacobianInverse ( 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::MappedGeometry< Map, Geo >::jacobian().

Referenced by Dune::MappedGeometry< Map, Geo >::jacobianInverseTransposed().

◆ jacobianInverseTransposed()

template<class Map , class Geo >
JacobianInverseTransposed Dune::MappedGeometry< Map, Geo >::jacobianInverseTransposed ( 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. \]

References Dune::MappedGeometry< Map, Geo >::jacobianInverse(), and Dune::transpose().

◆ jacobianTransposed()

template<class Map , class Geo >
JacobianTransposed Dune::MappedGeometry< Map, Geo >::jacobianTransposed ( const LocalCoordinate local) const
inline

Obtain the transposed of the Jacobian.

Parameters
[in]locallocal coordinate to evaluate Jacobian in
Returns
the matrix corresponding to the transposed of the Jacobian

References Dune::MappedGeometry< Map, Geo >::jacobian(), and Dune::transpose().

Referenced by Dune::MappedGeometry< Map, Geo >::integrationElement().

◆ LocalCoordinate()

template<class Map , class Geo >
Dune::MappedGeometry< Map, Geo >::LocalCoordinate ( const GlobalCoordinate y,
Impl::GaussNewtonOptions< ctype opts = {} 
) const
inline

Evaluate the inverse coordinate mapping.

Parameters
[in]yGlobal coordinate to map
[in]optsParameters to control the behavior of the Gauss-Newton algorithm.
Returns
For given global coordinate y the local coordinate x that minimizes the function (global(x) - y).two_norm2() over the local coordinate space spanned by the reference element.
Exceptions
Incase of an error indicating that the local coordinate can not be obtained, a Dune::Exception is thrown, with an error code from `GaussNewtonErrorCode`.
Note
It is not guaranteed that the resulting local coordinate is inside the reference element domain.

◆ volume()

template<class Map , class Geo >
Volume Dune::MappedGeometry< Map, Geo >::volume ( Impl::ConvergenceOptions< ctype opts = {}) const
inline

Obtain the volume of the mapping's image.

Calculates the volume of the entity by numerical integration. Since the polynomial order of the volume element is not known, iteratively compute numerical integrals with increasing order of the quadrature rules, until tolerance is reached.

Parameters
optsAn optional control over the convergence, providing a break tolerance and a maximal iteration count.

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.111.3 (Nov 21, 23:30, 2024)