|
dune-geometry 2.10
|
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. | |
| 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. | |
| LocalCoordinate | local (const GlobalCoordinate &y, Impl::GaussNewtonOptions< ctype > opts={}) const |
| Evaluate the inverse coordinate mapping. | |
| ctype | integrationElement (const LocalCoordinate &local) const |
| Obtain the integration element. | |
| Volume | volume (Impl::ConvergenceOptions< ctype > opts={}) const |
| Obtain the volume of the mapping's image. | |
| 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. | |
| JacobianTransposed | jacobianTransposed (const LocalCoordinate &local) const |
| Obtain the transposed of the Jacobian. | |
| JacobianInverse | jacobianInverse (const LocalCoordinate &local) const |
| Obtain the Jacobian's inverse. | |
| JacobianInverseTransposed | jacobianInverseTransposed (const LocalCoordinate &local) const |
| Obtain the transposed of the Jacobian's inverse. | |
Static Public Attributes | |
| static constexpr int | mydimension = LocalCoordinate::size() |
| geometry dimension | |
| static constexpr int | coorddimension = GlobalCoordinate::size() |
| coordinate dimension | |
Protected Types | |
| using | MatrixHelper = Impl::FieldMatrixHelper<ctype> |
| using | Mapping = Map |
| using | Geometry = Geo |
| using | DerivativeMapping = std::remove_reference_t<decltype(derivative(std::declval<Map>()))> |
Protected Member Functions | |
| ReferenceElement | refElement () const |
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.
The following expressions must be valid:
| Map | Differentiable mapping f: LF -> GF with LF = GG the range of the geometry g. |
| Geo | A geometry type g: LG -> GG fulfilling the concept Concept::Geometry. |
| using Dune::MappedGeometry< Map, Geo >::ctype = typename Geo::ctype |
coordinate type
|
protected |
|
protected |
| using Dune::MappedGeometry< Map, Geo >::GlobalCoordinate = std::remove_reference_t<decltype(std::declval<Map>()(std::declval<typename Geo::GlobalCoordinate>()))> |
type of global coordinates
| using Dune::MappedGeometry< Map, Geo >::Jacobian = FieldMatrix<ctype, coorddimension, mydimension> |
type of jacobian
| using Dune::MappedGeometry< Map, Geo >::JacobianInverse = FieldMatrix<ctype, mydimension, coorddimension> |
type of jacobian inverse
| using Dune::MappedGeometry< Map, Geo >::JacobianInverseTransposed = FieldMatrix<ctype, coorddimension, mydimension> |
type of jacobian inverse transposed
| using Dune::MappedGeometry< Map, Geo >::JacobianTransposed = FieldMatrix<ctype, mydimension, coorddimension> |
type of jacobian transposed
| using Dune::MappedGeometry< Map, Geo >::LocalCoordinate = typename Geo::LocalCoordinate |
type of local coordinates
|
protected |
|
protected |
| using Dune::MappedGeometry< Map, Geo >::Volume = std::remove_reference_t<decltype(Dune::power(std::declval<ctype>(),mydimension))> |
type of volume
|
inline |
Constructor from mapping to parametrize the geometry.
| [in] | mapping | A differentiable function for the parametrization of the geometry |
| [in] | geometry | The geometry that is mapped |
| [in] | affine | Flag indicating whether the MappedGeometry represents an affine mapping |
|
inline |
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.
|
inline |
Map the center of the wrapped geometry.
|
inline |
Obtain coordinates of the i-th corner.
|
inline |
Obtain number of corners of the corresponding reference element.
|
inline |
Evaluate the coordinate mapping.
Chained evaluation of the geometry and mapping from local coordinates in the reference element associated to the wrapped geometry.
| [in] | local | local coordinates |
|
inline |
Obtain the integration element.
If the Jacobian of the geometry is denoted by 

![\[ \mu(x) = \sqrt{|\det (J^T(x) J(x))|}.\]](form_1.png)
| [in] | local | local coordinate to evaluate the integration element in. |

|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Evaluate the inverse coordinate mapping.
| [in] | y | Global coordinate to map |
| [in] | opts | Parameters to control the behavior of the Gauss-Newton algorithm. |
| In | case of an error indicating that the local coordinate can not be obtained, a Dune::Exception is thrown, with an error code from `GaussNewtonErrorCode`. |
|
inlineprotected |
|
inline |
Obtain the geometry type from the reference element.
|
inline |
Obtain the volume of the mapping's image by given quadrature rules.
|
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.
| opts | An optional control over the convergence, providing a break tolerance and a maximal iteration count. |
|
staticconstexpr |
coordinate dimension
|
staticconstexpr |
geometry dimension