PTEMS  0.1.0-dev+git.81fd0e4
PolyTopic Element Method Solver
ptems::ReferenceMapping< DIM, SpaceDIM > Class Template Reference

Structure for defining a mapping between an element and a reference element. More...

#include <ptems/referencemapping.hpp>

Collaboration diagram for ptems::ReferenceMapping< DIM, SpaceDIM >:

Public Member Functions

constexpr ReferenceMapping ()
 Create an empty reference mapping. More...
 
constexpr ReferenceElement< DIM > Reference () const
 Gets the type of all extrusions from the reference element. More...
 
template<std::size_t D>
constexpr ExtrudeType Extrusion () const
 Gets the type of extrusion of the reference element in the \((D+1)\)th-dimension. More...
 
constexpr ExtrudeType Extrusion (std::size_t dim) const
 Gets the type of extrusion of the reference element in the \((dim+1)\)th-dimension. More...
 
Vector< SpaceDIM > LocalToGlobal (const Vector< DIM > &pt) const
 Maps a point in local coordinate space on the reference element to the global coordinate space. More...
 
template<std::size_t D = SpaceDIM, typename = std::enable_if_t<(D == SpaceDIM) && (D == DIM)>>
Vector< DIM > GlobalToLocal (const Vector< D > &pt) const
 Maps a point in global coordinate space to a point in the local coordinate space on the reference element. More...
 
Matrix< SpaceDIM, DIM > Jacobian (const Vector< DIM > &pt) const
 Gets the Jacobian matrix for the specified point. More...
 
std::array< Matrix< SpaceDIM >, DIM > Hessian (const Vector< DIM > &pt) const
 Gets the Hessian matrix for each component of the mapping for the specified point. More...
 
constexpr std::size_t VertexCount () const
 Gets the number of vertices in the reference element. More...
 
bool IsAffine () const
 Gets if the mapping is GUARANTEED to be affine. More...
 

Static Public Member Functions

template<ExtrudeType... Types>
constexpr static ReferenceMapping Create (const std::initializer_list< Vector< SpaceDIM >> &args)
 Creates mapping to the reference element defined by Types from the specified shape. More...
 

Static Public Attributes

static constexpr ReferenceElement< DIM > ReferenceSimplex = CreateReferenceElement(ExtrudeType::Point)
 Simplex reference element of dimension DIM. More...
 
static constexpr ReferenceElement< DIM > ReferenceHypercube = CreateReferenceElement(ExtrudeType::Hyperplane)
 Hypercube reference element of dimension DIM. More...
 

Detailed Description

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
class ptems::ReferenceMapping< DIM, SpaceDIM >

Structure for defining a mapping between an element and a reference element.

A 2D reference element can be constructed by considering the reference interval \([-1,1]\) in 1D and then "extruding" this between y=-1 and y=1 either as an interval at \(y=-1\) and \(y=1\) (the reference square \([-1,1]^2\)), or as an interval at \(y=-1\) and a point \((-1,1)\) at \(y=1\) (the reference triangle \((-1,-1)--(1,-1)--(-1,1)\)).

A 3D reference element can be constructed by extruding one of the two 2D reference elements at \(z=-1\) to either the same 2D reference element at \(z=1\), or to the point \((-1,-1,1)\):

  • Square extruded to square => reference cube \([-1,1]^3\)
  • Square extruded to point => reference pyramid \((-1,-1,-1),(1,-1,-1),(1,1,-1),(-1,1,-1),(-1,-1,1)\)
  • Triangle extruded to triangle => reference prism \((-1,-1,-1),(1,-1,-1),(-1,1,-1),(-1,-1,1),(1,-1,1),(-1,1,1)\)
  • Triangle extruded to point => reference tetrahedron \((-1,-1,-1),(1,-1,-1),(-1,1,-1),(-1,-1,1)\)

This can be repeated for any dimension. In general, therefore, we call the extrusion of the DIM-1 reference element to a DIM-1 reference element an extrusion to a "hyperplane" (ExtrudeType::Hyperplane), and to a point as extrusion to a "point" (ExtrudeType::Point)

Template Parameters
DIMDimension of the element being mapped
SpaceDIMThe dimension of the space (DIM <= SpaceD) the element is within

Constructor & Destructor Documentation

◆ ReferenceMapping()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr ptems::ReferenceMapping< DIM, SpaceDIM >::ReferenceMapping ( )
inlineconstexpr

Create an empty reference mapping.

Member Function Documentation

◆ Create()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
template<ExtrudeType... Types>
constexpr static ReferenceMapping ptems::ReferenceMapping< DIM, SpaceDIM >::Create ( const std::initializer_list< Vector< SpaceDIM >> &  args)
inlinestaticconstexpr

Creates mapping to the reference element defined by Types from the specified shape.

Template Parameters
TypesDIM-1 list of the extrusions for the \(N=2,3,...,DIM\) dimension to get the \(N\)-dimensional reference element from the \((N-1)\)-dimensional reference element, where the \(1\)-dimensional reference element is the interval \([-1,1]\).
Parameters
argsList of vertices of the element to map to the reference element. If the last extrusion in Types is ExtrudeType::Hyperplane, then the first half of vertices are the vertices on the \(\hat{x}_{DIM}=-1\) plane and the second half are the vertices on the \(\hat{x}_{DIM}=1\). If the last extrusion in Types is ExtrudeType::Point, then the last point is the single point on the \(\hat{x}_{DIM}=1\) plane and the rest are the vertices of the hyperplane on the \(\hat{x}_{DIM}=-1\) plane. (This is applied recursive for the \(N=DIM-1,...,1\) planes to each set of vertices defining a hyperplane)
Returns
The mapping from the specified element to the reference element

◆ Extrusion() [1/2]

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
template<std::size_t D>
constexpr ExtrudeType ptems::ReferenceMapping< DIM, SpaceDIM >::Extrusion ( ) const
inlineconstexpr

Gets the type of extrusion of the reference element in the \((D+1)\)th-dimension.

Template Parameters
DThe dimension to get the extrusion for
Returns
The extrusion type. Always ExtrudeType::Point if D=0.

◆ Extrusion() [2/2]

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr ExtrudeType ptems::ReferenceMapping< DIM, SpaceDIM >::Extrusion ( std::size_t  dim) const
inlineconstexpr

Gets the type of extrusion of the reference element in the \((dim+1)\)th-dimension.

Exceptions
std::out_of_rangeIf dim is not less than DIM
Template Parameters
dimThe dimension to get the extrusion for
Returns
The extrusion type. Always ExtrudeType::Point if dim=0.

◆ GlobalToLocal()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
template<std::size_t D = SpaceDIM, typename = std::enable_if_t<(D == SpaceDIM) && (D == DIM)>>
Vector<DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::GlobalToLocal ( const Vector< D > &  pt) const
inline

Maps a point in global coordinate space to a point in the local coordinate space on the reference element.

Only exists if SpaceDIM == DIM

Parameters
ptPoint in global coordinate to map
Returns
The point mapped into local coordinates space on the reference element

◆ Hessian()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
std::array<Matrix<SpaceDIM>, DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::Hessian ( const Vector< DIM > &  pt) const
inline

Gets the Hessian matrix for each component of the mapping for the specified point.

Parameters
ptPoint in local coordinate space to get Hessian for
Returns
The Hessian for each component of the mapping

◆ IsAffine()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
bool ptems::ReferenceMapping< DIM, SpaceDIM >::IsAffine ( ) const
inline

Gets if the mapping is GUARANTEED to be affine.

Warning
This method may return false even if the mapping is affine. If this function returns true then the mapping is guaranteed to be affine. Essentially,

\[ \textrm{IsAffine()} = \textrm{true} \quad \implies \quad \text{Mapping Affine} \]

but

\[ \text{Mapping Affine} \quad \not\implies \quad \textrm{IsAffine()} = \textrm{true} \]

Note
The current implementation returns true if and only if the mapping is for the reference simplex.
Returns
true If the mapping is definitely affine; false if it may be non-affine.

◆ Jacobian()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
Matrix<SpaceDIM, DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::Jacobian ( const Vector< DIM > &  pt) const
inline

Gets the Jacobian matrix for the specified point.

Note
If IsAffine() is true then this method will return the same matrix for ANY point.
Parameters
ptPoint in local coordinate space to get Jacobian for
Returns
The Jacobian for the specified point

◆ LocalToGlobal()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
Vector<SpaceDIM> ptems::ReferenceMapping< DIM, SpaceDIM >::LocalToGlobal ( const Vector< DIM > &  pt) const
inline

Maps a point in local coordinate space on the reference element to the global coordinate space.

Parameters
ptPoint in local coordinate space on the reference element to map
Returns
The point mapped into global coordinates

◆ Reference()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr ReferenceElement<DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::Reference ( ) const
inlineconstexpr

Gets the type of all extrusions from the reference element.

Returns
The extrusion types defining the reference element

◆ VertexCount()

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr std::size_t ptems::ReferenceMapping< DIM, SpaceDIM >::VertexCount ( ) const
inlineconstexpr

Gets the number of vertices in the reference element.

Returns
Number of vertices in reference element

Member Data Documentation

◆ ReferenceHypercube

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr ReferenceElement<DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::ReferenceHypercube = CreateReferenceElement(ExtrudeType::Hyperplane)
staticconstexpr

Hypercube reference element of dimension DIM.

◆ ReferenceSimplex

template<std::size_t DIM, std::size_t SpaceDIM = DIM>
constexpr ReferenceElement<DIM> ptems::ReferenceMapping< DIM, SpaceDIM >::ReferenceSimplex = CreateReferenceElement(ExtrudeType::Point)
staticconstexpr

Simplex reference element of dimension DIM.


The documentation for this class was generated from the following file: