MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::DofTransformation Class Referenceabstract

#include <doftrans.hpp>

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

Public Member Functions

int Size () const
 
int Height () const
 
int NumRows () const
 
int Width () const
 
int NumCols () const
 
void SetFaceOrientations (const Array< int > &face_orientation)
 Configure the transformation using face orientations for the current element. More...
 
const Array< int > & GetFaceOrientations () const
 
virtual void TransformPrimal (double *v) const =0
 
virtual void TransformPrimal (Vector &v) const
 
virtual void TransformPrimalCols (DenseMatrix &V) const
 Transform groups of DoFs stored as dense matrices. More...
 
virtual void InvTransformPrimal (double *v) const =0
 
virtual void InvTransformPrimal (Vector &v) const
 
virtual void TransformDual (double *v) const =0
 
virtual void TransformDual (Vector &v) const
 
virtual void InvTransformDual (double *v) const =0
 
virtual void InvTransformDual (Vector &v) const
 
virtual void TransformDual (DenseMatrix &V) const
 
virtual void TransformDualRows (DenseMatrix &V) const
 Transform groups of dual DoFs stored as dense matrices. More...
 
virtual void TransformDualCols (DenseMatrix &V) const
 
virtual ~DofTransformation ()
 

Protected Member Functions

 DofTransformation (int size)
 

Protected Attributes

int size_
 
Array< int > Fo
 

Detailed Description

The DofTransformation 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. These transformations are necessary to ensure that basis functions in neighboring elements align correctly. Closely related but complementary transformations are required for the entries stored in LinearForm and BilinearForm objects. The DofTransformation 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 56 of file doftrans.hpp.

Constructor & Destructor Documentation

mfem::DofTransformation::DofTransformation ( int  size)
inlineprotected

Definition at line 63 of file doftrans.hpp.

virtual mfem::DofTransformation::~DofTransformation ( )
inlinevirtual

Definition at line 117 of file doftrans.hpp.

Member Function Documentation

const Array<int>& mfem::DofTransformation::GetFaceOrientations ( ) const
inline

Definition at line 80 of file doftrans.hpp.

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

Definition at line 69 of file doftrans.hpp.

virtual void mfem::DofTransformation::InvTransformDual ( double *  v) const
pure virtual
void mfem::DofTransformation::InvTransformDual ( Vector v) const
virtual

Definition at line 88 of file doftrans.cpp.

virtual void mfem::DofTransformation::InvTransformPrimal ( 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_WedgeDofTransformation, mfem::ND_TetDofTransformation, mfem::ND_TriDofTransformation, and mfem::VDofTransformation.

void mfem::DofTransformation::InvTransformPrimal ( Vector v) const
virtual

Definition at line 60 of file doftrans.cpp.

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

Definition at line 72 of file doftrans.hpp.

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

Definition at line 70 of file doftrans.hpp.

void mfem::DofTransformation::SetFaceOrientations ( const Array< int > &  face_orientation)
inline

Configure the transformation using face orientations for the current element.

The face_orientation array can be obtained from Mesh::GetElementFaces.

Definition at line 77 of file doftrans.hpp.

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

Definition at line 68 of file doftrans.hpp.

virtual void mfem::DofTransformation::TransformDual ( double *  v) const
pure virtual

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

Implemented in mfem::ND_WedgeDofTransformation, mfem::ND_TetDofTransformation, mfem::ND_TriDofTransformation, and mfem::VDofTransformation.

void mfem::DofTransformation::TransformDual ( Vector v) const
virtual

Definition at line 30 of file doftrans.cpp.

void mfem::DofTransformation::TransformDual ( DenseMatrix V) const
virtual

Transform a matrix of dual DoFs entries as computed by a BilinearFormIntegrator before summing into a BilinearForm object.

Definition at line 35 of file doftrans.cpp.

void mfem::DofTransformation::TransformDualCols ( DenseMatrix V) const
virtual

Definition at line 52 of file doftrans.cpp.

void mfem::DofTransformation::TransformDualRows ( DenseMatrix V) const
virtual

Transform groups of dual DoFs stored as dense matrices.

Definition at line 41 of file doftrans.cpp.

virtual void mfem::DofTransformation::TransformPrimal ( 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_WedgeDofTransformation, mfem::ND_TetDofTransformation, mfem::ND_TriDofTransformation, and mfem::VDofTransformation.

void mfem::DofTransformation::TransformPrimal ( Vector v) const
virtual

Definition at line 17 of file doftrans.cpp.

void mfem::DofTransformation::TransformPrimalCols ( DenseMatrix V) const
virtual

Transform groups of DoFs stored as dense matrices.

Definition at line 22 of file doftrans.cpp.

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

Definition at line 71 of file doftrans.hpp.

Member Data Documentation

Array<int> mfem::DofTransformation::Fo
protected

Definition at line 61 of file doftrans.hpp.

int mfem::DofTransformation::size_
protected

Definition at line 59 of file doftrans.hpp.


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