MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::IsoparametricTransformation Class Reference

A standard isoparametric element transformation. More...

#include <eltrans.hpp>

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

Public Member Functions

 IsoparametricTransformation ()
 
void SetFE (const FiniteElement *FE)
 Set the element that will be used to compute the transformations.
 
const FiniteElementGetFE () const
 Get the current element used to compute the transformations.
 
void SetPointMat (const DenseMatrix &pm)
 Set the underlying point matrix describing the transformation.
 
const DenseMatrixGetPointMat () const
 Return the stored point matrix.
 
DenseMatrixGetPointMat ()
 Write access to the stored point matrix. Use with caution.
 
void SetIdentityTransformation (Geometry::Type GeomType)
 Set the FiniteElement Geometry for the reference elements being used.
 
virtual void Transform (const IntegrationPoint &, Vector &)
 Transform integration point from reference coordinates to physical coordinates and store them in the vector.
 
virtual void Transform (const IntegrationRule &, DenseMatrix &)
 Transform all the integration points from the integration rule from reference coordinates to physical coordinates and store them as column vectors in the matrix.
 
virtual void Transform (const DenseMatrix &matrix, DenseMatrix &result)
 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.
 
virtual int Order () const
 Return the order of the current element we are using for the transformation.
 
virtual int OrderJ () const
 Return the order of the elements of the Jacobian of the transformation.
 
virtual int OrderW () const
 Return the order of the determinant of the Jacobian (weight) of the transformation.
 
virtual int OrderGrad (const FiniteElement *fe) const
 Return the order of \( adj(J)^T \nabla fi \).
 
virtual int GetSpaceDim () const
 Get the dimension of the target (physical) space.
 
virtual int TransformBack (const Vector &v, IntegrationPoint &ip, const real_t phys_rel_tol=1e-15)
 Transform a point pt from physical space to a point ip in reference space and optionally can set a solver tolerance using phys_tol.
 
virtual ~IsoparametricTransformation ()
 
MFEM_DEPRECATED void FinalizeTransformation ()
 
- Public Member Functions inherited from mfem::ElementTransformation
 ElementTransformation ()
 
void Reset ()
 Force the reevaluation of the Jacobian in the next call.
 
void SetIntPoint (const IntegrationPoint *ip)
 Set the integration point ip that weights and Jacobians will be evaluated at.
 
const IntegrationPointGetIntPoint ()
 Get a const reference to the currently set integration point. This will return NULL if no integration point is set.
 
const DenseMatrixJacobian ()
 Return the Jacobian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
 
const DenseMatrixHessian ()
 Return the Hessian matrix of the transformation at the currently set IntegrationPoint, using the method SetIntPoint().
 
real_t 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} \).
 
const DenseMatrixAdjugateJacobian ()
 Return the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint.
 
const DenseMatrixTransposeAdjugateJacobian ()
 Return the transpose of the adjugate of the Jacobian matrix of the transformation at the currently set IntegrationPoint.
 
const DenseMatrixInverseJacobian ()
 Return the inverse of the Jacobian matrix of the transformation at the currently set IntegrationPoint.
 
Geometry::Type GetGeometryType () const
 Return the Geometry::Type of the reference element.
 
int GetDimension () const
 Return the topological dimension of the reference element.
 
virtual ~ElementTransformation ()
 

Additional Inherited Members

- Public Types inherited from mfem::ElementTransformation
enum  {
  ELEMENT = 1 , BDR_ELEMENT = 2 , EDGE = 3 , FACE = 4 ,
  BDR_FACE = 5
}
 
- Public Attributes inherited from mfem::ElementTransformation
int Attribute
 
int ElementNo
 
int ElementType
 
const Meshmesh
 The Mesh object containing the element.
 
- Protected Types inherited from mfem::ElementTransformation
enum  StateMasks {
  JACOBIAN_MASK = 1 , WEIGHT_MASK = 2 , ADJUGATE_MASK = 4 , INVERSE_MASK = 8 ,
  HESSIAN_MASK = 16 , TRANS_ADJUGATE_MASK = 32
}
 
- Protected Member Functions inherited from mfem::ElementTransformation
real_t EvalWeight ()
 
const DenseMatrixEvalAdjugateJ ()
 
const DenseMatrixEvalTransAdjugateJ ()
 
const DenseMatrixEvalInverseJ ()
 
- Protected Attributes inherited from mfem::ElementTransformation
const IntegrationPointIntPoint
 
DenseMatrix dFdx
 
DenseMatrix adjJ
 
DenseMatrix invJ
 
DenseMatrix d2Fdx2
 
DenseMatrix adjJT
 
real_t Wght
 
int EvalState
 
Geometry::Type geom
 

Detailed Description

A standard isoparametric element transformation.

Definition at line 362 of file eltrans.hpp.

Constructor & Destructor Documentation

◆ IsoparametricTransformation()

mfem::IsoparametricTransformation::IsoparametricTransformation ( )
inline

Definition at line 379 of file eltrans.hpp.

◆ ~IsoparametricTransformation()

virtual mfem::IsoparametricTransformation::~IsoparametricTransformation ( )
inlinevirtual

Definition at line 459 of file eltrans.hpp.

Member Function Documentation

◆ FinalizeTransformation()

MFEM_DEPRECATED void mfem::IsoparametricTransformation::FinalizeTransformation ( )
inline

Definition at line 461 of file eltrans.hpp.

◆ GetFE()

const FiniteElement * mfem::IsoparametricTransformation::GetFE ( ) const
inline

Get the current element used to compute the transformations.

Definition at line 390 of file eltrans.hpp.

◆ GetPointMat() [1/2]

DenseMatrix & mfem::IsoparametricTransformation::GetPointMat ( )
inline

Write access to the stored point matrix. Use with caution.

If the point matrix is altered using this member function the Reset function should also be called to force the reevaluation of the Jacobian, etc..

Definition at line 411 of file eltrans.hpp.

◆ GetPointMat() [2/2]

const DenseMatrix & mfem::IsoparametricTransformation::GetPointMat ( ) const
inline

Return the stored point matrix.

Definition at line 405 of file eltrans.hpp.

◆ GetSpaceDim()

virtual int mfem::IsoparametricTransformation::GetSpaceDim ( ) const
inlinevirtual

Get the dimension of the target (physical) space.

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

Implements mfem::ElementTransformation.

Definition at line 443 of file eltrans.hpp.

◆ Order()

virtual int mfem::IsoparametricTransformation::Order ( ) const
inlinevirtual

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

Implements mfem::ElementTransformation.

Definition at line 431 of file eltrans.hpp.

◆ OrderGrad()

int mfem::IsoparametricTransformation::OrderGrad ( const FiniteElement * fe) const
virtual

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

Implements mfem::ElementTransformation.

Definition at line 467 of file eltrans.cpp.

◆ OrderJ()

int mfem::IsoparametricTransformation::OrderJ ( ) const
virtual

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

Implements mfem::ElementTransformation.

Definition at line 439 of file eltrans.cpp.

◆ OrderW()

int mfem::IsoparametricTransformation::OrderW ( ) const
virtual

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

Implements mfem::ElementTransformation.

Definition at line 453 of file eltrans.cpp.

◆ SetFE()

void mfem::IsoparametricTransformation::SetFE ( const FiniteElement * FE)
inline

Set the element that will be used to compute the transformations.

Definition at line 382 of file eltrans.hpp.

◆ SetIdentityTransformation()

void mfem::IsoparametricTransformation::SetIdentityTransformation ( Geometry::Type GeomType)

Set the FiniteElement Geometry for the reference elements being used.

Definition at line 379 of file eltrans.cpp.

◆ SetPointMat()

void mfem::IsoparametricTransformation::SetPointMat ( const DenseMatrix & pm)
inline

Set the underlying point matrix describing the transformation.

The dimensions of the matrix are space-dim x dof. The transformation is defined as \( x = F( \hat x ) = P \phi( \hat x ) \)

where \( \hat x \) is the reference point, x is the corresponding physical point, P is the point matrix, and \( \phi( \hat x ) \) is the column-vector of all basis functions evaluated at \( \hat x \) . The columns of P represent the control points in physical space defining the transformation.

Definition at line 402 of file eltrans.hpp.

◆ Transform() [1/3]

void mfem::IsoparametricTransformation::Transform ( const DenseMatrix & matrix,
DenseMatrix & result )
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.

Implements mfem::ElementTransformation.

Reimplemented in mfem::FaceElementTransformations.

Definition at line 525 of file eltrans.cpp.

◆ Transform() [2/3]

void mfem::IsoparametricTransformation::Transform ( const IntegrationPoint & ip,
Vector & trans )
virtual

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

Implements mfem::ElementTransformation.

Reimplemented in mfem::FaceElementTransformations.

Definition at line 488 of file eltrans.cpp.

◆ Transform() [3/3]

void mfem::IsoparametricTransformation::Transform ( const IntegrationRule & ir,
DenseMatrix & tr )
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.

Implements mfem::ElementTransformation.

Reimplemented in mfem::FaceElementTransformations.

Definition at line 499 of file eltrans.cpp.

◆ TransformBack()

virtual int mfem::IsoparametricTransformation::TransformBack ( const Vector & v,
IntegrationPoint & ip,
const real_t phys_rel_tol = 1e-15 )
inlinevirtual

Transform a point pt from physical space to a point ip in reference space and optionally can set a solver tolerance using phys_tol.

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.

Implements mfem::ElementTransformation.

Definition at line 451 of file eltrans.hpp.


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