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

#include <doftrans.hpp>

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

Public Member Functions

bool IsIdentity () const override
 If the DofTransformation performs no transformation.
 
void TransformPrimal (const Array< int > &Fo, real_t *v) const override
 
void InvTransformPrimal (const Array< int > &Fo, real_t *v) const override
 
void TransformDual (const Array< int > &Fo, real_t *v) const override
 
void InvTransformDual (const Array< int > &Fo, real_t *v) const override
 
- 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
 
virtual ~StatelessDofTransformation ()=default
 

Static Public Member Functions

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

Protected Member Functions

 ND_DofTransformation (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 301 of file doftrans.hpp.

Constructor & Destructor Documentation

◆ ND_DofTransformation()

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

Definition at line 203 of file doftrans.cpp.

Member Function Documentation

◆ GetFaceInverseTransform()

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

Definition at line 322 of file doftrans.hpp.

◆ GetFaceTransform()

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

Definition at line 319 of file doftrans.hpp.

◆ InvTransformDual()

void mfem::ND_DofTransformation::InvTransformDual ( const Array< int > & face_orientation,
real_t * v ) const
overridevirtual

Inverse Transform dual DoFs

Implements mfem::StatelessDofTransformation.

Definition at line 291 of file doftrans.cpp.

◆ InvTransformPrimal()

void mfem::ND_DofTransformation::InvTransformPrimal ( const Array< int > & face_orientation,
real_t * v ) const
overridevirtual

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 240 of file doftrans.cpp.

◆ IsIdentity()

bool mfem::ND_DofTransformation::IsIdentity ( ) const
inlineoverridevirtual

If the DofTransformation performs no transformation.

Implements mfem::StatelessDofTransformation.

Definition at line 325 of file doftrans.hpp.

◆ TransformDual()

void mfem::ND_DofTransformation::TransformDual ( const Array< int > & face_orientation,
real_t * v ) const
overridevirtual

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

Implements mfem::StatelessDofTransformation.

Definition at line 266 of file doftrans.cpp.

◆ TransformPrimal()

void mfem::ND_DofTransformation::TransformPrimal ( const Array< int > & face_orientation,
real_t * v ) const
overridevirtual

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 214 of file doftrans.cpp.

Member Data Documentation

◆ nedges

const int mfem::ND_DofTransformation::nedges
protected

Definition at line 312 of file doftrans.hpp.

◆ nedofs

const int mfem::ND_DofTransformation::nedofs
protected

Definition at line 310 of file doftrans.hpp.

◆ nfaces

const int mfem::ND_DofTransformation::nfaces
protected

Definition at line 313 of file doftrans.hpp.

◆ nfdofs

const int mfem::ND_DofTransformation::nfdofs
protected

Definition at line 311 of file doftrans.hpp.

◆ order

const int mfem::ND_DofTransformation::order
protected

Definition at line 309 of file doftrans.hpp.


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