MFEM  v4.6.0
Finite element discretization library
Public Types | Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | List of all members
mfem::ElementTransformation Class Referenceabstract

#include <eltrans.hpp>

Inheritance diagram for mfem::ElementTransformation:
[legend]
Collaboration diagram for mfem::ElementTransformation:
[legend]

Public Types

enum  {
  ELEMENT = 1, BDR_ELEMENT = 2, EDGE = 3, FACE = 4,
  BDR_FACE = 5
}
 

Public Member Functions

 ElementTransformation ()
 
void Reset ()
 Force the reevaluation of the Jacobian in the next call. More...
 
void SetIntPoint (const IntegrationPoint *ip)
 Set the integration point ip that weights and Jacobians will be evaluated at. More...
 
const IntegrationPointGetIntPoint ()
 Get a const reference to the currently set integration point. This will return NULL if no integration point is set. More...
 
virtual void Transform (const IntegrationPoint &, Vector &)=0
 Transform integration point from reference coordinates to physical coordinates and store them in the vector. More...
 
virtual void Transform (const IntegrationRule &, DenseMatrix &)=0
 Transform all the integration points from the integration rule from reference coordinates to physical coordinates and store them as column vectors in the matrix. More...
 
virtual void Transform (const DenseMatrix &matrix, DenseMatrix &result)=0
 Transform all the integration points from the column vectors of matrix from reference coordinates to physical coordinates and store them as column vectors in result. More...
 
const DenseMatrixJacobian ()
 Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint(). More...
 
const DenseMatrixHessian ()
 Return the Hessian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint(). More...
 
double Weight ()
 Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint. The Weight evaluates to \( \sqrt{\lvert J^T J \rvert} \). More...
 
const DenseMatrixAdjugateJacobian ()
 Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint. More...
 
const DenseMatrixTransposeAdjugateJacobian ()
 Return the transpose of the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint. More...
 
const DenseMatrixInverseJacobian ()
 Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint. More...
 
virtual int Order () const =0
 Return the order of the current element we are using for the transformation. More...
 
virtual int OrderJ () const =0
 Return the order of the elements of the Jacobian of the transformation. More...
 
virtual int OrderW () const =0
 Return the order of the determinant of the Jacobian (weight) of the transformation. More...
 
virtual int OrderGrad (const FiniteElement *fe) const =0
 Return the order of \( adj(J)^T \nabla fi \). More...
 
Geometry::Type GetGeometryType () const
 Return the Geometry::Type of the reference element. More...
 
int GetDimension () const
 Return the topological dimension of the reference element. More...
 
virtual int GetSpaceDim () const =0
 Get the dimension of the target (physical) space. More...
 
virtual int TransformBack (const Vector &pt, IntegrationPoint &ip)=0
 Transform a point pt from physical space to a point ip in reference space. More...
 
virtual ~ElementTransformation ()
 

Public Attributes

int Attribute
 
int ElementNo
 
int ElementType
 
class Meshmesh
 The Mesh object containing the element. More...
 

Protected Types

enum  StateMasks {
  JACOBIAN_MASK = 1, WEIGHT_MASK = 2, ADJUGATE_MASK = 4, INVERSE_MASK = 8,
  HESSIAN_MASK = 16, TRANS_ADJUGATE_MASK = 32
}
 

Protected Member Functions

virtual const DenseMatrixEvalJacobian ()=0
 Evaluate the Jacobian of the transformation at the IntPoint and store it in dFdx. More...
 
virtual const DenseMatrixEvalHessian ()=0
 Evaluate the Hessian of the transformation at the IntPoint and store it in d2Fdx2. More...
 
double EvalWeight ()
 
const DenseMatrixEvalAdjugateJ ()
 
const DenseMatrixEvalTransAdjugateJ ()
 
const DenseMatrixEvalInverseJ ()
 

Protected Attributes

const IntegrationPointIntPoint
 
DenseMatrix dFdx
 
DenseMatrix adjJ
 
DenseMatrix invJ
 
DenseMatrix d2Fdx2
 
DenseMatrix adjJT
 
double Wght
 
int EvalState
 
Geometry::Type geom
 

Detailed Description

Definition at line 23 of file eltrans.hpp.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum

This enumeration declares the values stored in ElementTransformation::ElementType and indicates which group of objects the index stored in ElementTransformation::ElementNo refers:

| ElementType | Range of ElementNo +----------—+----------------------— | ELEMENT | [0, Mesh::GetNE() ) | BDR_ELEMENT | [0, Mesh::GetNBE() ) | EDGE | [0, Mesh::GetNEdges() ) | FACE | [0, Mesh::GetNFaces() ) | BDR_FACE | [0, Mesh::GetNBE() )

Enumerator
ELEMENT 
BDR_ELEMENT 
EDGE 
FACE 
BDR_FACE 

Definition at line 69 of file eltrans.hpp.

◆ StateMasks

Enumerator
JACOBIAN_MASK 
WEIGHT_MASK 
ADJUGATE_MASK 
INVERSE_MASK 
HESSIAN_MASK 
TRANS_ADJUGATE_MASK 

Definition at line 31 of file eltrans.hpp.

Constructor & Destructor Documentation

◆ ElementTransformation()

mfem::ElementTransformation::ElementTransformation ( )

Definition at line 19 of file eltrans.cpp.

◆ ~ElementTransformation()

virtual mfem::ElementTransformation::~ElementTransformation ( )
inlinevirtual

Definition at line 180 of file eltrans.hpp.

Member Function Documentation

◆ AdjugateJacobian()

const DenseMatrix& mfem::ElementTransformation::AdjugateJacobian ( )
inline

Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint.

Definition at line 135 of file eltrans.hpp.

◆ EvalAdjugateJ()

const DenseMatrix & mfem::ElementTransformation::EvalAdjugateJ ( )
protected

Definition at line 36 of file eltrans.cpp.

◆ EvalHessian()

virtual const DenseMatrix& mfem::ElementTransformation::EvalHessian ( )
protectedpure virtual

Evaluate the Hessian of the transformation at the IntPoint and store it in d2Fdx2.

◆ EvalInverseJ()

const DenseMatrix & mfem::ElementTransformation::EvalInverseJ ( )
protected

Definition at line 57 of file eltrans.cpp.

◆ EvalJacobian()

virtual const DenseMatrix& mfem::ElementTransformation::EvalJacobian ( )
protectedpure virtual

Evaluate the Jacobian of the transformation at the IntPoint and store it in dFdx.

◆ EvalTransAdjugateJ()

const DenseMatrix & mfem::ElementTransformation::EvalTransAdjugateJ ( )
protected

Definition at line 46 of file eltrans.cpp.

◆ EvalWeight()

double mfem::ElementTransformation::EvalWeight ( )
protected

Definition at line 28 of file eltrans.cpp.

◆ GetDimension()

int mfem::ElementTransformation::GetDimension ( ) const
inline

Return the topological dimension of the reference element.

Definition at line 165 of file eltrans.hpp.

◆ GetGeometryType()

Geometry::Type mfem::ElementTransformation::GetGeometryType ( ) const
inline

Return the Geometry::Type of the reference element.

Definition at line 162 of file eltrans.hpp.

◆ GetIntPoint()

const IntegrationPoint& mfem::ElementTransformation::GetIntPoint ( )
inline

Get a const reference to the currently set integration point. This will return NULL if no integration point is set.

Definition at line 98 of file eltrans.hpp.

◆ GetSpaceDim()

virtual int mfem::ElementTransformation::GetSpaceDim ( ) const
pure virtual

Get the dimension of the target (physical) space.

We support 2D meshes embedded in 3D; in this case the function will return "3".

Implemented in mfem::IsoparametricTransformation.

◆ Hessian()

const DenseMatrix& mfem::ElementTransformation::Hessian ( )
inline

Return the Hessian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().

Definition at line 125 of file eltrans.hpp.

◆ InverseJacobian()

const DenseMatrix& mfem::ElementTransformation::InverseJacobian ( )
inline

Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint.

Definition at line 145 of file eltrans.hpp.

◆ Jacobian()

const DenseMatrix& mfem::ElementTransformation::Jacobian ( )
inline

Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().

The dimensions of the Jacobian matrix are physical-space-dim by reference-space-dim. The first column contains the x derivatives of the transformation, the second – the y derivatives, etc.

Definition at line 119 of file eltrans.hpp.

◆ Order()

virtual int mfem::ElementTransformation::Order ( ) const
pure virtual

Return the order of the current element we are using for the transformation.

Implemented in mfem::IsoparametricTransformation.

◆ OrderGrad()

virtual int mfem::ElementTransformation::OrderGrad ( const FiniteElement fe) const
pure virtual

Return the order of \( adj(J)^T \nabla fi \).

Implemented in mfem::IsoparametricTransformation.

◆ OrderJ()

virtual int mfem::ElementTransformation::OrderJ ( ) const
pure virtual

Return the order of the elements of the Jacobian of the transformation.

Implemented in mfem::IsoparametricTransformation.

◆ OrderW()

virtual int mfem::ElementTransformation::OrderW ( ) const
pure virtual

Return the order of the determinant of the Jacobian (weight) of the transformation.

Implemented in mfem::IsoparametricTransformation.

◆ Reset()

void mfem::ElementTransformation::Reset ( )
inline

Force the reevaluation of the Jacobian in the next call.

Definition at line 89 of file eltrans.hpp.

◆ SetIntPoint()

void mfem::ElementTransformation::SetIntPoint ( const IntegrationPoint ip)
inline

Set the integration point ip that weights and Jacobians will be evaluated at.

Definition at line 93 of file eltrans.hpp.

◆ Transform() [1/3]

virtual void mfem::ElementTransformation::Transform ( const IntegrationPoint ,
Vector  
)
pure virtual

Transform integration point from reference coordinates to physical coordinates and store them in the vector.

Implemented in mfem::FaceElementTransformations, and mfem::IsoparametricTransformation.

◆ Transform() [2/3]

virtual void mfem::ElementTransformation::Transform ( const IntegrationRule ,
DenseMatrix  
)
pure virtual

Transform all the integration points from the integration rule from reference coordinates to physical coordinates and store them as column vectors in the matrix.

Implemented in mfem::FaceElementTransformations, and mfem::IsoparametricTransformation.

◆ Transform() [3/3]

virtual void mfem::ElementTransformation::Transform ( const DenseMatrix matrix,
DenseMatrix result 
)
pure virtual

Transform all the integration points from the column vectors of matrix from reference coordinates to physical coordinates and store them as column vectors in result.

Implemented in mfem::FaceElementTransformations, and mfem::IsoparametricTransformation.

◆ TransformBack()

virtual int mfem::ElementTransformation::TransformBack ( const Vector pt,
IntegrationPoint ip 
)
pure virtual

Transform a point pt from physical space to a point ip in reference space.

Attempt to find the IntegrationPoint that is transformed into the given point in physical space. If the inversion fails a non-zero value is returned. This method is not 100 percent reliable for non-linear transformations.

Implemented in mfem::IsoparametricTransformation.

◆ TransposeAdjugateJacobian()

const DenseMatrix& mfem::ElementTransformation::TransposeAdjugateJacobian ( )
inline

Return the transpose of the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint.

Definition at line 140 of file eltrans.hpp.

◆ Weight()

double mfem::ElementTransformation::Weight ( )
inline

Return the weight of the Jacobian matrix of the transformation at the currently set IntegrationPoint. The Weight evaluates to \( \sqrt{\lvert J^T J \rvert} \).

Definition at line 131 of file eltrans.hpp.

Member Data Documentation

◆ adjJ

DenseMatrix mfem::ElementTransformation::adjJ
protected

Definition at line 27 of file eltrans.hpp.

◆ adjJT

DenseMatrix mfem::ElementTransformation::adjJT
protected

Definition at line 28 of file eltrans.hpp.

◆ Attribute

int mfem::ElementTransformation::Attribute

Definition at line 78 of file eltrans.hpp.

◆ d2Fdx2

DenseMatrix mfem::ElementTransformation::d2Fdx2
protected

Definition at line 28 of file eltrans.hpp.

◆ dFdx

DenseMatrix mfem::ElementTransformation::dFdx
protected

Definition at line 27 of file eltrans.hpp.

◆ ElementNo

int mfem::ElementTransformation::ElementNo

Definition at line 78 of file eltrans.hpp.

◆ ElementType

int mfem::ElementTransformation::ElementType

Definition at line 78 of file eltrans.hpp.

◆ EvalState

int mfem::ElementTransformation::EvalState
protected

Definition at line 30 of file eltrans.hpp.

◆ geom

Geometry::Type mfem::ElementTransformation::geom
protected

Definition at line 40 of file eltrans.hpp.

◆ IntPoint

const IntegrationPoint* mfem::ElementTransformation::IntPoint
protected

Definition at line 26 of file eltrans.hpp.

◆ invJ

DenseMatrix mfem::ElementTransformation::invJ
protected

Definition at line 27 of file eltrans.hpp.

◆ mesh

class Mesh* mfem::ElementTransformation::mesh

The Mesh object containing the element.

If the element transformation belongs to a mesh, this will point to the containing Mesh object. ElementNo will be the number of the element in this Mesh. This will be NULL if the element does not belong to a mesh.

Definition at line 84 of file eltrans.hpp.

◆ Wght

double mfem::ElementTransformation::Wght
protected

Definition at line 29 of file eltrans.hpp.


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