MFEM v4.8.0
Finite element discretization library
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
mfem::FaceElementTransformations Class Reference

A specialized ElementTransformation class representing a face and its two neighboring elements. More...

#include <eltrans.hpp>

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

Public Types

enum  ConfigMasks {
  HAVE_ELEM1 = 1 , HAVE_ELEM2 = 2 , HAVE_LOC1 = 4 , HAVE_LOC2 = 8 ,
  HAVE_FACE = 16
}
 
- Public Types inherited from mfem::ElementTransformation
enum  {
  ELEMENT = 1 , BDR_ELEMENT = 2 , EDGE = 3 , FACE = 4 ,
  BDR_FACE = 5
}
 

Public Member Functions

 FaceElementTransformations ()
 
void SetGeometryType (Geometry::Type g)
 Method to set the geometry type of the face.
 
int GetConfigurationMask () const
 Return the mask defining the configuration state.
 
void SetIntPoint (const IntegrationPoint *face_ip)
 Set the integration point in the Face and the two neighboring elements, if present.
 
void SetAllIntPoints (const IntegrationPoint *face_ip)
 Set the integration point in the Face and the two neighboring elements, if present.
 
const IntegrationPointGetElement1IntPoint ()
 Get a const reference to the integration point in neighboring element 1 corresponding to the currently set integration point on the face.
 
const IntegrationPointGetElement2IntPoint ()
 Get a const reference to the integration point in neighboring element 2 corresponding to the currently set integration point on the face.
 
void Transform (const IntegrationPoint &, Vector &) override
 Transform integration point from reference coordinates to physical coordinates and store them in the vector.
 
void Transform (const IntegrationRule &, DenseMatrix &) override
 Transform all the integration points from the integration rule from reference coordinates to physical coordinates and store them as column vectors in the matrix.
 
void Transform (const DenseMatrix &matrix, DenseMatrix &result) override
 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.
 
ElementTransformationGetElement1Transformation ()
 
ElementTransformationGetElement2Transformation ()
 
IntegrationPointTransformationGetIntPoint1Transformation ()
 
IntegrationPointTransformationGetIntPoint2Transformation ()
 
real_t CheckConsistency (int print_level=0, std::ostream &out=mfem::out)
 Check for self-consistency: compares the result of mapping the reference face vertices to physical coordinates using the three transformations: face, element 1, and element 2.
 
- Public Member Functions inherited from mfem::IsoparametricTransformation
 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.
 
void Transform (const IntegrationPoint &, Vector &) override
 Transform integration point from reference coordinates to physical coordinates and store them in the vector.
 
void Transform (const IntegrationRule &, DenseMatrix &) override
 Transform all the integration points from the integration rule from reference coordinates to physical coordinates and store them as column vectors in the matrix.
 
void Transform (const DenseMatrix &matrix, DenseMatrix &result) override
 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.
 
int Order () const override
 Return the order of the current element we are using for the transformation.
 
int OrderJ () const override
 Return the order of the elements of the Jacobian of the transformation.
 
int OrderW () const override
 Return the order of the determinant of the Jacobian (weight) of the transformation.
 
int OrderGrad (const FiniteElement *fe) const override
 Return the order of adj(J)Tfi.
 
int GetSpaceDim () const override
 Get the dimension of the target (physical) space.
 
int TransformBack (const Vector &v, IntegrationPoint &ip, const real_t phys_rel_tol=tol_0) override
 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 |JTJ|.
 
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 ()
 

Public Attributes

int Elem1No
 
int Elem2No
 
Geometry::TypeFaceGeom
 
ElementTransformationElem1
 
ElementTransformationElem2
 
ElementTransformationFace
 
IntegrationPointTransformation Loc1
 
IntegrationPointTransformation Loc2
 
- Public Attributes inherited from mfem::ElementTransformation
int Attribute
 
int ElementNo
 
int ElementType
 
const Meshmesh
 The Mesh object containing the element.
 

Protected Member Functions

void SetConfigurationMask (int m)
 Set the mask indicating which portions of the object have been setup.
 
- Protected Member Functions inherited from mfem::ElementTransformation
real_t EvalWeight ()
 
const DenseMatrixEvalAdjugateJ ()
 
const DenseMatrixEvalTransAdjugateJ ()
 
const DenseMatrixEvalInverseJ ()
 

Friends

class Mesh
 
class ParMesh
 

Additional Inherited Members

- 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 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
 
- Static Protected Attributes inherited from mfem::ElementTransformation
static constexpr real_t tol_0 = 1e-15
 

Detailed Description

A specialized ElementTransformation class representing a face and its two neighboring elements.

This class can be used as a container for the element transformation data needed for integrating discontinuous fields on element interfaces in a Discontinuous Galerkin (DG) context.

The secondary purpose of this class is to enable the GridFunction::GetValue function, and various related functions, to properly evaluate fields with limited continuity on boundary elements.

Definition at line 749 of file eltrans.hpp.

Member Enumeration Documentation

◆ ConfigMasks

Enumerator
HAVE_ELEM1 

Element on side 1 is configured.

HAVE_ELEM2 

Element on side 2 is configured.

HAVE_LOC1 

Point transformation for side 1 is configured.

HAVE_LOC2 

Point transformation for side 2 is configured.

HAVE_FACE 

Face transformation is configured.

Definition at line 780 of file eltrans.hpp.

Constructor & Destructor Documentation

◆ FaceElementTransformations()

mfem::FaceElementTransformations::FaceElementTransformations ( )
inline

Definition at line 795 of file eltrans.hpp.

Member Function Documentation

◆ CheckConsistency()

real_t mfem::FaceElementTransformations::CheckConsistency ( int print_level = 0,
std::ostream & out = mfem::out )

Check for self-consistency: compares the result of mapping the reference face vertices to physical coordinates using the three transformations: face, element 1, and element 2.

Parameters
[in]print_levelIf set to a positive number, print the physical coordinates of the face vertices computed through all available transformations: face, element 1, and/or element 2.
[in,out]outThe output stream to use for printing.
Returns
A maximal distance between physical coordinates of face vertices that should coincide. A successful check should return a small number relative to the mesh extents. If less than 2 of the three transformations are set, returns 0.
Warning
This check will generally fail on periodic boundary faces.

Definition at line 687 of file eltrans.cpp.

◆ GetConfigurationMask()

int mfem::FaceElementTransformations::GetConfigurationMask ( ) const
inline

Return the mask defining the configuration state.

The mask value indicates which portions of FaceElementTransformations object have been configured.

mask & 1: Elem1 is configured mask & 2: Elem2 is configured mask & 4: Loc1 is configured mask & 8: Loc2 is configured mask & 16: The Face transformation itself is configured

Definition at line 818 of file eltrans.hpp.

◆ GetElement1IntPoint()

const IntegrationPoint & mfem::FaceElementTransformations::GetElement1IntPoint ( )
inline

Get a const reference to the integration point in neighboring element 1 corresponding to the currently set integration point on the face.

This IntegrationPoint object will only contain up-to-date data if SetIntPoint or SetAllIntPoints has been called with the latest integration point for the face and the appropriate point transformation has been configured.

Definition at line 846 of file eltrans.hpp.

◆ GetElement1Transformation()

ElementTransformation & mfem::FaceElementTransformations::GetElement1Transformation ( )

Definition at line 632 of file eltrans.cpp.

◆ GetElement2IntPoint()

const IntegrationPoint & mfem::FaceElementTransformations::GetElement2IntPoint ( )
inline

Get a const reference to the integration point in neighboring element 2 corresponding to the currently set integration point on the face.

This IntegrationPoint object will only contain up-to-date data if SetIntPoint or SetAllIntPoints has been called with the latest integration point for the face and the appropriate point transformation has been configured.

Definition at line 856 of file eltrans.hpp.

◆ GetElement2Transformation()

ElementTransformation & mfem::FaceElementTransformations::GetElement2Transformation ( )

Definition at line 640 of file eltrans.cpp.

◆ GetIntPoint1Transformation()

IntegrationPointTransformation & mfem::FaceElementTransformations::GetIntPoint1Transformation ( )

Definition at line 648 of file eltrans.cpp.

◆ GetIntPoint2Transformation()

IntegrationPointTransformation & mfem::FaceElementTransformations::GetIntPoint2Transformation ( )

Definition at line 656 of file eltrans.cpp.

◆ SetAllIntPoints()

void mfem::FaceElementTransformations::SetAllIntPoints ( const IntegrationPoint * face_ip)
inline

Set the integration point in the Face and the two neighboring elements, if present.

This is a more expressive member function name than SetIntPoint, which in this special case, does the same thing. This function can be used for greater code clarity.

Definition at line 835 of file eltrans.hpp.

◆ SetConfigurationMask()

void mfem::FaceElementTransformations::SetConfigurationMask ( int m)
inlineprotected

Set the mask indicating which portions of the object have been setup.

The argument m is a bitmask used in Mesh::GetFaceElementTransformations to indicate which portions of the FaceElementTransformations object have been configured.

mask & 1: Elem1 is configured mask & 2: Elem2 is configured mask & 4: Loc1 is configured mask & 8: Loc2 is configured mask & 16: The Face transformation itself is configured

Definition at line 776 of file eltrans.hpp.

◆ SetGeometryType()

void mfem::FaceElementTransformations::SetGeometryType ( Geometry::Type g)
inline

Method to set the geometry type of the face.

Note
This method is designed to be used when [Par]Mesh::GetFaceTransformation will not be called i.e. when the face transformation will not be needed but the neighboring element transformations will be. Using this method to override the GeometryType should only be done with great care.

Definition at line 805 of file eltrans.hpp.

◆ SetIntPoint()

void mfem::FaceElementTransformations::SetIntPoint ( const IntegrationPoint * face_ip)

Set the integration point in the Face and the two neighboring elements, if present.

The point face_ip must be in the reference coordinate system of the face.

Definition at line 609 of file eltrans.cpp.

◆ Transform() [1/3]

void mfem::FaceElementTransformations::Transform ( const DenseMatrix & matrix,
DenseMatrix & result )
overridevirtual

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.

Definition at line 679 of file eltrans.cpp.

◆ Transform() [2/3]

void mfem::FaceElementTransformations::Transform ( const IntegrationPoint & ,
Vector &  )
overridevirtual

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

Implements mfem::ElementTransformation.

Definition at line 663 of file eltrans.cpp.

◆ Transform() [3/3]

void mfem::FaceElementTransformations::Transform ( const IntegrationRule & ,
DenseMatrix &  )
overridevirtual

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.

Definition at line 671 of file eltrans.cpp.

Friends And Related Symbol Documentation

◆ Mesh

friend class Mesh
friend

Definition at line 760 of file eltrans.hpp.

◆ ParMesh

friend class ParMesh
friend

Definition at line 762 of file eltrans.hpp.

Member Data Documentation

◆ Elem1

ElementTransformation* mfem::FaceElementTransformations::Elem1

Definition at line 791 of file eltrans.hpp.

◆ Elem1No

int mfem::FaceElementTransformations::Elem1No

Definition at line 789 of file eltrans.hpp.

◆ Elem2

ElementTransformation * mfem::FaceElementTransformations::Elem2

Definition at line 791 of file eltrans.hpp.

◆ Elem2No

int mfem::FaceElementTransformations::Elem2No

Definition at line 789 of file eltrans.hpp.

◆ Face

ElementTransformation* mfem::FaceElementTransformations::Face
Deprecated
No longer necessary

Definition at line 792 of file eltrans.hpp.

◆ FaceGeom

Geometry::Type& mfem::FaceElementTransformations::FaceGeom
Deprecated
Use GetGeometryType instead

Definition at line 790 of file eltrans.hpp.

◆ Loc1

IntegrationPointTransformation mfem::FaceElementTransformations::Loc1

Definition at line 793 of file eltrans.hpp.

◆ Loc2

IntegrationPointTransformation mfem::FaceElementTransformations::Loc2

Definition at line 793 of file eltrans.hpp.


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