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

Transfer data between a coarse mesh and an embedded refined mesh using interpolation. More...

#include <fespace.hpp>

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

Public Member Functions

 InterpolationGridTransfer (FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
 
virtual ~InterpolationGridTransfer ()
 
void SetMassIntegrator (BilinearFormIntegrator *mass_integ_, bool own_mass_integ_=true)
 Assign a mass integrator to be used in the construction of the backward, fine-to-coarse, transfer operator. More...
 
virtual const OperatorForwardOperator ()
 Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the range FE space. More...
 
virtual const OperatorBackwardOperator ()
 Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the domain FE space. More...
 
- Public Member Functions inherited from mfem::GridTransfer
 GridTransfer (FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
 
virtual ~GridTransfer ()
 Virtual destructor. More...
 
void SetOperatorType (Operator::Type type)
 Set the desired Operator::Type for the construction of all operators defined by the underlying transfer algorithm. More...
 
virtual const OperatorTrueForwardOperator ()
 Return an Operator that transfers true-dof Vectors from the domain FE space to true-dof Vectors in the range FE space. More...
 
virtual const OperatorTrueBackwardOperator ()
 Return an Operator that transfers true-dof Vectors from the range FE space back to true-dof Vectors in the domain FE space. More...
 

Protected Attributes

BilinearFormIntegratormass_integ
 Ownership depends on own_mass_integ. More...
 
bool own_mass_integ
 Ownership flag for mass_integ. More...
 
OperatorHandle F
 Forward, coarse-to-fine, operator. More...
 
OperatorHandle B
 Backward, fine-to-coarse, operator. More...
 
- Protected Attributes inherited from mfem::GridTransfer
FiniteElementSpacedom_fes
 Domain FE space. More...
 
FiniteElementSpaceran_fes
 Range FE space. More...
 
Operator::Type oper_type
 Desired Operator::Type for the construction of all operators defined by the underlying transfer algorithm. It can be ignored by derived classes. More...
 
OperatorHandle fw_t_oper
 Forward true-dof operator. More...
 
OperatorHandle bw_t_oper
 Backward true-dof operator. More...
 
bool parallel
 

Additional Inherited Members

- Protected Member Functions inherited from mfem::GridTransfer
bool Parallel () const
 
const OperatorMakeTrueOperator (FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
 

Detailed Description

Transfer data between a coarse mesh and an embedded refined mesh using interpolation.

The forward, coarse-to-fine, transfer uses nodal interpolation. The backward, fine-to-coarse, transfer is defined locally (on a coarse element) as B = (F^t M_f F)^{-1} F^t M_f, where F is the forward transfer matrix, and M_f is a mass matrix on the union of all fine elements comprising the coarse element. Note that the backward transfer operator, B, is a left inverse of the forward transfer operator, F, i.e. B F = I. Both F and B are defined in reference space and do not depend on the actual physical shape of the mesh elements.

It is assumed that both the coarse and the fine FiniteElementSpaces use compatible types of elements, e.g. finite elements with the same map-type (VALUE, INTEGRAL, H_DIV, H_CURL - see class FiniteElement). Generally, the FE spaces can have different orders, however, in order for the backward operator to be well-defined, the (local) number of the fine dofs should not be smaller than the number of coarse dofs.

Definition at line 807 of file fespace.hpp.

Constructor & Destructor Documentation

mfem::InterpolationGridTransfer::InterpolationGridTransfer ( FiniteElementSpace coarse_fes,
FiniteElementSpace fine_fes 
)
inline

Definition at line 817 of file fespace.hpp.

mfem::InterpolationGridTransfer::~InterpolationGridTransfer ( )
virtual

Definition at line 2480 of file fespace.cpp.

Member Function Documentation

const Operator & mfem::InterpolationGridTransfer::BackwardOperator ( )
virtual

Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the domain FE space.

Implements mfem::GridTransfer.

Definition at line 2527 of file fespace.cpp.

const Operator & mfem::InterpolationGridTransfer::ForwardOperator ( )
virtual

Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the range FE space.

Implements mfem::GridTransfer.

Definition at line 2494 of file fespace.cpp.

void mfem::InterpolationGridTransfer::SetMassIntegrator ( BilinearFormIntegrator mass_integ_,
bool  own_mass_integ_ = true 
)

Assign a mass integrator to be used in the construction of the backward, fine-to-coarse, transfer operator.

Definition at line 2485 of file fespace.cpp.

Member Data Documentation

OperatorHandle mfem::InterpolationGridTransfer::B
protected

Backward, fine-to-coarse, operator.

Definition at line 814 of file fespace.hpp.

OperatorHandle mfem::InterpolationGridTransfer::F
protected

Forward, coarse-to-fine, operator.

Definition at line 813 of file fespace.hpp.

BilinearFormIntegrator* mfem::InterpolationGridTransfer::mass_integ
protected

Ownership depends on own_mass_integ.

Definition at line 810 of file fespace.hpp.

bool mfem::InterpolationGridTransfer::own_mass_integ
protected

Ownership flag for mass_integ.

Definition at line 811 of file fespace.hpp.


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