Dune Core Modules (unstable)

Dune::LocalFiniteElementGeometry< LFE, cdim > Class Template Reference

Geometry implementation based on local-basis function parametrization. More...

#include <dune/geometry/localfiniteelementgeometry.hh>

Public Types

using ctype = typename LocalBasisTraits::DomainFieldType
 coordinate type
 
using LocalCoordinate = FieldVector< ctype, mydimension >
 type of local coordinates
 
using GlobalCoordinate = FieldVector< ctype, coorddimension >
 type of global coordinates
 
using Volume = decltype(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
 
using ReferenceElements = Dune::ReferenceElements< ctype, mydimension >
 type of reference element
 

Public Member Functions

 LocalFiniteElementGeometry ()=default
 Default constructed geometry results in an empty/invalid representation.
 
 LocalFiniteElementGeometry (const ReferenceElement &refElement, const LocalFiniteElement &localFE, std::vector< GlobalCoordinate > vertices)
 Constructor from a vector of coefficients of the LocalBasis parametrizing the geometry. More...
 
template<class Param , std::enable_if_t< std::is_invocable_r_v< GlobalCoordinate, Param, LocalCoordinate >, int > = 0>
 LocalFiniteElementGeometry (const ReferenceElement &refElement, const LocalFiniteElement &localFE, Param &&parametrization)
 Constructor from a local parametrization function, mapping local to (curved) global coordinates. More...
 
template<class... Args>
 LocalFiniteElementGeometry (GeometryType gt, Args &&... args)
 Constructor, forwarding to the other constructors that take a reference-element. More...
 
int order () const
 Obtain the polynomial order of the parametrization.
 
bool affine () const
 Is this mapping affine? This is only true for flat affine geometries.
 
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 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...
 
const LocalFiniteElement & finiteElement () const
 Obtain the local finite-element.
 
const std::vector< GlobalCoordinate > & coefficients () const
 Obtain the coefficients of the parametrization.
 
const LocalBasis & localBasis () const
 The local basis of the stored local finite-element.
 

Static Public Attributes

static const int mydimension = LocalBasisTraits::dimDomain
 geometry dimension
 
static const int coorddimension = cdim
 coordinate dimension
 

Detailed Description

template<class LFE, int cdim>
class Dune::LocalFiniteElementGeometry< LFE, cdim >

Geometry implementation based on local-basis function parametrization.

Parametrization of the geometry by any localfunction interpolated into a local finite-element space.

Template Parameters
LFEType of a local finite-element.
cdimCoordinate dimension.

Constructor & Destructor Documentation

◆ LocalFiniteElementGeometry() [1/3]

template<class LFE , int cdim>
Dune::LocalFiniteElementGeometry< LFE, cdim >::LocalFiniteElementGeometry ( const ReferenceElement &  refElement,
const LocalFiniteElement &  localFE,
std::vector< GlobalCoordinate vertices 
)
inline

Constructor from a vector of coefficients of the LocalBasis parametrizing the geometry.

Parameters
[in]refElementreference element for the geometry
[in]localFELocal finite-element to use for the parametrization
[in]verticesCoefficients of the local interpolation into the basis

The vertices are the coefficients of a local interpolation in the local finite-element. For Lagrange local bases, these correspond to vertices on the curved geometry in the local Lagrange nodes.

Note
The vertices are stored internally, so if possible move an external vertex storage to this constructor

◆ LocalFiniteElementGeometry() [2/3]

template<class LFE , int cdim>
template<class Param , std::enable_if_t< std::is_invocable_r_v< GlobalCoordinate, Param, LocalCoordinate >, int > = 0>
Dune::LocalFiniteElementGeometry< LFE, cdim >::LocalFiniteElementGeometry ( const ReferenceElement &  refElement,
const LocalFiniteElement &  localFE,
Param &&  parametrization 
)
inline

Constructor from a local parametrization function, mapping local to (curved) global coordinates.

Parameters
[in]refElementreference element for the geometry
[in]localFELocal finite-element to use for the parametrization
[in]parametrizationparametrization function with signature GlobalCoordinate(LocalCoordinate)

The parametrization function is not stored in the class, but interpolated into the local finite-element basis and the computed interpolation coefficients are stored.

◆ LocalFiniteElementGeometry() [3/3]

template<class LFE , int cdim>
template<class... Args>
Dune::LocalFiniteElementGeometry< LFE, cdim >::LocalFiniteElementGeometry ( GeometryType  gt,
Args &&...  args 
)
inlineexplicit

Constructor, forwarding to the other constructors that take a reference-element.

Parameters
[in]gtgeometry type
[in]args...arguments passed to the other constructors

Member Function Documentation

◆ global()

template<class LFE , int cdim>
GlobalCoordinate Dune::LocalFiniteElementGeometry< LFE, cdim >::global ( const LocalCoordinate local) const
inline

Evaluate the coordinate mapping.

Implements a linear combination of local basis functions scaled by the vertices as coefficients.

\[ global = \sum_i v_i \psi_i(local) \]

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

References Dune::DenseVector< V >::axpy(), and Dune::LocalFiniteElementGeometry< LFE, cdim >::localBasis().

Referenced by Dune::LocalFiniteElementGeometry< LFE, cdim >::center(), and Dune::LocalFiniteElementGeometry< LFE, cdim >::corner().

◆ integrationElement()

template<class LFE , int cdim>
ctype Dune::LocalFiniteElementGeometry< LFE, cdim >::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::LocalFiniteElementGeometry< LFE, cdim >::jacobianTransposed().

Referenced by Dune::LocalFiniteElementGeometry< LFE, cdim >::volume().

◆ jacobian()

template<class LFE , int cdim>
Jacobian Dune::LocalFiniteElementGeometry< LFE, cdim >::jacobian ( const LocalCoordinate local) const
inline

Obtain the Jacobian.

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

References Dune::LocalFiniteElementGeometry< LFE, cdim >::localBasis(), and Dune::FieldMatrix< ctype, coorddimension, mydimension >::rows.

Referenced by Dune::LocalFiniteElementGeometry< LFE, cdim >::jacobianInverse(), and Dune::LocalFiniteElementGeometry< LFE, cdim >::jacobianTransposed().

◆ jacobianInverse()

template<class LFE , int cdim>
JacobianInverse Dune::LocalFiniteElementGeometry< LFE, cdim >::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::LocalFiniteElementGeometry< LFE, cdim >::jacobian().

Referenced by Dune::LocalFiniteElementGeometry< LFE, cdim >::jacobianInverseTransposed().

◆ jacobianInverseTransposed()

template<class LFE , int cdim>
JacobianInverseTransposed Dune::LocalFiniteElementGeometry< LFE, cdim >::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::LocalFiniteElementGeometry< LFE, cdim >::jacobianInverse(), and Dune::FieldMatrix< K, ROWS, COLS >::transposed().

◆ jacobianTransposed()

template<class LFE , int cdim>
JacobianTransposed Dune::LocalFiniteElementGeometry< LFE, cdim >::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::LocalFiniteElementGeometry< LFE, cdim >::jacobian(), and Dune::FieldMatrix< K, ROWS, COLS >::transposed().

Referenced by Dune::LocalFiniteElementGeometry< LFE, cdim >::integrationElement().

◆ LocalCoordinate()

template<class LFE , int cdim>
Dune::LocalFiniteElementGeometry< LFE, cdim >::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, an 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 LFE , int cdim>
Volume Dune::LocalFiniteElementGeometry< LFE, cdim >::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.


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 (Jul 15, 22:36, 2024)