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

#include <doftrans.hpp>

Inheritance diagram for mfem::StatelessDofTransformation:
[legend]

Public Member Functions

int Size () const
 
int Height () const
 
int NumRows () const
 
int Width () const
 
int NumCols () const
 
virtual void TransformPrimal (const Array< int > &face_orientation, double *v) const =0
 
void TransformPrimal (const Array< int > &face_orientation, Vector &v) const
 
virtual void InvTransformPrimal (const Array< int > &face_orientation, double *v) const =0
 
void InvTransformPrimal (const Array< int > &face_orientation, Vector &v) const
 
virtual void TransformDual (const Array< int > &face_orientation, double *v) const =0
 
void TransformDual (const Array< int > &face_orientation, Vector &v) const
 
virtual void InvTransformDual (const Array< int > &face_orientation, double *v) const =0
 
void InvTransformDual (const Array< int > &face_orientation, Vector &v) const
 

Protected Member Functions

 StatelessDofTransformation (int size)
 

Protected Attributes

int size_
 

Detailed Description

The StatelessDofTransformation class is an abstract base class for a family of transformations that map local degrees of freedom (DoFs), contained within individual elements, to global degrees of freedom, stored within GridFunction objects.

In this context "stateless" means that the concrete classes derived from StatelessDofTransformation do not store information about the relative orientations of the faces with respect to their neighboring elements. In other words there is no information specific to a particular element (aside from the element type e.g. tetrahedron, wedge, or pyramid). The StatelessDofTransformation provides access to the transformation operators for specific relative face orientations. These are useful, for example, when relating DoFs associated with distinct overlapping meshes such as parent and sub-meshes.

These transformations are necessary to ensure that basis functions in neighboring (or overlapping) elements align correctly. Closely related but complementary transformations are required for the entries stored in LinearForm and BilinearForm objects. The StatelessDofTransformation class is designed to apply the action of both of these types of DoF transformations.

Let the "primal transformation" be given by the operator T. This means that given a local element vector v the data that must be placed into a GridFunction object is v_t = T * v.

We also need the inverse of the primal transformation T^{-1} so that we can recover the local element vector from data read out of a GridFunction e.g. v = T^{-1} * v_t.

We need to preserve the action of our linear forms applied to primal vectors. In other words, if f is the local vector computed by a linear form then f * v = f_t * v_t (where "*" represents an inner product of vectors). This requires that f_t = T^{-T} * f i.e. the "dual transform" is given by the transpose of the inverse of the primal transformation.

For bilinear forms we require that v^T * A * v = v_t^T * A_t * v_t. This implies that A_t = T^{-T} * A * T^{-1}. This can be accomplished by performing dual transformations of the rows and columns of the matrix A.

For discrete linear operators the range must be modified with the primal transformation rather than the dual transformation because the result is a primal vector rather than a dual vector. This leads to the transformation D_t = T * D * T^{-1}. This can be accomplished by using a primal transformation on the columns of D and a dual transformation on its rows.

Definition at line 68 of file doftrans.hpp.

Constructor & Destructor Documentation

◆ StatelessDofTransformation()

mfem::StatelessDofTransformation::StatelessDofTransformation ( int  size)
inlineprotected

Definition at line 73 of file doftrans.hpp.

Member Function Documentation

◆ Height()

int mfem::StatelessDofTransformation::Height ( ) const
inline

Definition at line 78 of file doftrans.hpp.

◆ InvTransformDual() [1/2]

virtual void mfem::StatelessDofTransformation::InvTransformDual ( const Array< int > &  face_orientation,
double *  v 
) const
pure virtual

Inverse Transform dual DoFs

Implemented in mfem::ND_StatelessDofTransformation, and mfem::StatelessVDofTransformation.

◆ InvTransformDual() [2/2]

void mfem::StatelessDofTransformation::InvTransformDual ( const Array< int > &  face_orientation,
Vector v 
) const
inline

Definition at line 115 of file doftrans.hpp.

◆ InvTransformPrimal() [1/2]

virtual void mfem::StatelessDofTransformation::InvTransformPrimal ( const Array< int > &  face_orientation,
double *  v 
) const
pure virtual

Inverse transform local DoFs. Used to transform DoFs from a global vector back to their element-local form. For example, this must be used to transform the vector obtained using GridFunction::GetSubVector before it can be used to compute a local interpolation.

Implemented in mfem::ND_StatelessDofTransformation, and mfem::StatelessVDofTransformation.

◆ InvTransformPrimal() [2/2]

void mfem::StatelessDofTransformation::InvTransformPrimal ( const Array< int > &  face_orientation,
Vector v 
) const
inline

Definition at line 100 of file doftrans.hpp.

◆ NumCols()

int mfem::StatelessDofTransformation::NumCols ( ) const
inline

Definition at line 81 of file doftrans.hpp.

◆ NumRows()

int mfem::StatelessDofTransformation::NumRows ( ) const
inline

Definition at line 79 of file doftrans.hpp.

◆ Size()

int mfem::StatelessDofTransformation::Size ( ) const
inline

Definition at line 77 of file doftrans.hpp.

◆ TransformDual() [1/2]

virtual void mfem::StatelessDofTransformation::TransformDual ( const Array< int > &  face_orientation,
double *  v 
) const
pure virtual

Transform dual DoFs as computed by a LinearFormIntegrator before summing into a LinearForm object.

Implemented in mfem::ND_StatelessDofTransformation, and mfem::StatelessVDofTransformation.

◆ TransformDual() [2/2]

void mfem::StatelessDofTransformation::TransformDual ( const Array< int > &  face_orientation,
Vector v 
) const
inline

Definition at line 108 of file doftrans.hpp.

◆ TransformPrimal() [1/2]

virtual void mfem::StatelessDofTransformation::TransformPrimal ( const Array< int > &  face_orientation,
double *  v 
) const
pure virtual

Transform local DoFs to align with the global DoFs. For example, this transformation can be used to map the local vector computed by FiniteElement::Project() to the transformed vector stored within a GridFunction object.

Implemented in mfem::ND_StatelessDofTransformation, and mfem::StatelessVDofTransformation.

◆ TransformPrimal() [2/2]

void mfem::StatelessDofTransformation::TransformPrimal ( const Array< int > &  face_orientation,
Vector v 
) const
inline

Definition at line 89 of file doftrans.hpp.

◆ Width()

int mfem::StatelessDofTransformation::Width ( ) const
inline

Definition at line 80 of file doftrans.hpp.

Member Data Documentation

◆ size_

int mfem::StatelessDofTransformation::size_
protected

Definition at line 71 of file doftrans.hpp.


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