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

#include <doftrans.hpp>

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

Public Member Functions

void TransformPrimal (const Array< int > &face_orientation, double *v) const
 
void InvTransformPrimal (const Array< int > &face_orientation, double *v) const
 
void TransformDual (const Array< int > &face_orientation, double *v) const
 
void InvTransformDual (const Array< int > &face_orientation, double *v) const
 
- Public Member Functions inherited from mfem::StatelessDofTransformation
int Size () const
 
int Height () const
 
int NumRows () const
 
int Width () const
 
int NumCols () const
 
void TransformPrimal (const Array< int > &face_orientation, Vector &v) const
 
void InvTransformPrimal (const Array< int > &face_orientation, Vector &v) const
 
void TransformDual (const Array< int > &face_orientation, Vector &v) const
 
void InvTransformDual (const Array< int > &face_orientation, Vector &v) const
 

Static Public Member Functions

static const DenseMatrixGetFaceTransform (int ori)
 
static const DenseMatrixGetFaceInverseTransform (int ori)
 

Protected Member Functions

 ND_StatelessDofTransformation (int size, int order, int num_edges, int num_tri_faces)
 
- Protected Member Functions inherited from mfem::StatelessDofTransformation
 StatelessDofTransformation (int size)
 

Protected Attributes

const int order
 
const int nedofs
 
const int nfdofs
 
const int nedges
 
const int nfaces
 
- Protected Attributes inherited from mfem::StatelessDofTransformation
int size_
 

Detailed Description

Abstract base class for high-order Nedelec spaces on elements with triangular faces.

The Nedelec DoFs on the interior of triangular faces come in pairs which share an interpolation point but have different vector directions. These directions depend on the orientation of the face and can therefore differ in neighboring elements. The mapping required to transform these DoFs can be implemented as series of 2x2 linear transformations. The raw data for these linear transformations is stored in the T_data and TInv_data arrays and can be accessed as DenseMatrices using the GetFaceTransform() and GetFaceInverseTransform() methods.

Definition at line 399 of file doftrans.hpp.

Constructor & Destructor Documentation

◆ ND_StatelessDofTransformation()

mfem::ND_StatelessDofTransformation::ND_StatelessDofTransformation ( int  size,
int  order,
int  num_edges,
int  num_tri_faces 
)
protected

Definition at line 212 of file doftrans.cpp.

Member Function Documentation

◆ GetFaceInverseTransform()

static const DenseMatrix& mfem::ND_StatelessDofTransformation::GetFaceInverseTransform ( int  ori)
inlinestatic

Definition at line 421 of file doftrans.hpp.

◆ GetFaceTransform()

static const DenseMatrix& mfem::ND_StatelessDofTransformation::GetFaceTransform ( int  ori)
inlinestatic

Definition at line 418 of file doftrans.hpp.

◆ InvTransformDual()

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

Inverse Transform dual DoFs

Implements mfem::StatelessDofTransformation.

Definition at line 296 of file doftrans.cpp.

◆ InvTransformPrimal()

void mfem::ND_StatelessDofTransformation::InvTransformPrimal ( const Array< int > &  face_orientation,
double *  v 
) const
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.

Implements mfem::StatelessDofTransformation.

Definition at line 248 of file doftrans.cpp.

◆ TransformDual()

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

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

Implements mfem::StatelessDofTransformation.

Definition at line 272 of file doftrans.cpp.

◆ TransformPrimal()

void mfem::ND_StatelessDofTransformation::TransformPrimal ( const Array< int > &  face_orientation,
double *  v 
) const
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.

Implements mfem::StatelessDofTransformation.

Definition at line 224 of file doftrans.cpp.

Member Data Documentation

◆ nedges

const int mfem::ND_StatelessDofTransformation::nedges
protected

Definition at line 410 of file doftrans.hpp.

◆ nedofs

const int mfem::ND_StatelessDofTransformation::nedofs
protected

Definition at line 408 of file doftrans.hpp.

◆ nfaces

const int mfem::ND_StatelessDofTransformation::nfaces
protected

Definition at line 411 of file doftrans.hpp.

◆ nfdofs

const int mfem::ND_StatelessDofTransformation::nfdofs
protected

Definition at line 409 of file doftrans.hpp.

◆ order

const int mfem::ND_StatelessDofTransformation::order
protected

Definition at line 407 of file doftrans.hpp.


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