MFEM
v4.2.0
Finite element discretization library
|
Transfer data between a coarse mesh and an embedded refined mesh using interpolation. More...
#include <fespace.hpp>
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 Operator & | ForwardOperator () |
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the range FE space. More... | |
virtual const Operator & | BackwardOperator () |
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 Operator & | TrueForwardOperator () |
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 Operator & | TrueBackwardOperator () |
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 | |
BilinearFormIntegrator * | mass_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 | |
FiniteElementSpace & | dom_fes |
Domain FE space. More... | |
FiniteElementSpace & | ran_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 Operator & | MakeTrueOperator (FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper) |
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 870 of file fespace.hpp.
|
inline |
Definition at line 880 of file fespace.hpp.
|
virtual |
Definition at line 2656 of file fespace.cpp.
|
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 2703 of file fespace.cpp.
|
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 2670 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 2661 of file fespace.cpp.
|
protected |
Backward, fine-to-coarse, operator.
Definition at line 877 of file fespace.hpp.
|
protected |
Forward, coarse-to-fine, operator.
Definition at line 876 of file fespace.hpp.
|
protected |
Ownership depends on own_mass_integ.
Definition at line 873 of file fespace.hpp.
|
protected |
Ownership flag for mass_integ.
Definition at line 874 of file fespace.hpp.