12 #ifndef MFEM_TRANSFER_HPP 13 #define MFEM_TRANSFER_HPP 15 #include "../linalg/linalg.hpp" 144 bool own_mass_integ_ =
true);
185 virtual void SetRelTol(
double p_rtol_) = 0;
186 virtual void SetAbsTol(
double p_atol_) = 0;
342 bool force_l2_space_ =
false)
431 const Operator* elem_restrict_lex_l;
432 const Operator* elem_restrict_lex_h;
Abstract class for all finite elements.
virtual ~InterpolationGridTransfer()
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
virtual void Prolongate(const Vector &x, Vector &y) const
Conjugate gradient method.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Data type for scaled Jacobi-type smoother of sparse matrix.
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
L2ProjectionH1Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
L2Prolongation(const L2Projection &l2proj_)
virtual void Mult(const Vector &x, Vector &y) const
FiniteElementSpace & dom_fes
Domain FE space.
Matrix-free transfer operator between finite element spaces on the same mesh.
Pointer to an Operator of a specified type.
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
virtual ~L2Prolongation()
Data type dense matrix using column-major storage.
virtual ~TransferOperator()
Destructor.
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
virtual void ProlongateTranspose(const Vector &x, Vector &y) const =0
~TrueTransferOperator()
Destructor.
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
virtual void Prolongate(const Vector &x, Vector &y) const
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
void BuildHo2Lor(int nel_ho, int nel_lor, const CoarseFineTransformations &cf_tr)
virtual ~PRefinementTransferOperator()
Destructor.
virtual const Operator & TrueBackwardOperator()
Return an Operator that transfers true-dof Vectors from the range FE space back to true-dof Vectors i...
virtual ~GridTransfer()
Virtual destructor.
FiniteElementSpace & ran_fes
Range FE space.
const FiniteElementSpace & fes_ho
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
virtual void SetRelTol(double p_rtol_)=0
void ElemMixedMass(Geometry::Type geom, const FiniteElement &fe_ho, const FiniteElement &fe_lor, ElementTransformation *el_tr, IntegrationPointTransformation &ip_tr, DenseMatrix &M_mixed_el) const
OperatorHandle F
Forward, coarse-to-fine, operator.
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
virtual void Prolongate(const Vector &x, Vector &y) const =0
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual void SetAbsTol(double p_atol_)
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
OperatorHandle bw_t_oper
Backward true-dof operator.
virtual bool SupportsBackwardsOperator() const
virtual void Mult(const Vector &x, Vector &y) const
void SetOperatorType(Operator::Type type)
Set the desired Operator::Type for the construction of all operators defined by the underlying transf...
OperatorHandle B
Backward, fine-to-coarse, operator.
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
OperatorHandle fw_t_oper
Forward true-dof operator.
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.
Type
Enumeration defining IDs for some classes derived from Operator.
void MultTranspose(const Vector &x, Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
Matrix-free transfer operator between finite element spaces.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
Matrix-free transfer operator between finite element spaces on the same mesh exploiting the tensor pr...
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
L2Prolongation * B
Backward, fine-to-coarse, operator.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
bool own_mass_integ
Ownership flag for mass_integ.
virtual void MultTranspose(const Vector &x, Vector &y) const
L2Projection * F
Forward, coarse-to-fine, operator.
virtual void SetAbsTol(double p_atol_)=0
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
virtual const Operator & TrueForwardOperator()
Return an Operator that transfers true-dof Vectors from the domain FE space to true-dof Vectors in th...
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a vector x corresponding to the coarse space to the vector y corresp...
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes_, FiniteElementSpace &fine_fes_, bool force_l2_space_=false)
const FiniteElementSpace & fes_lor
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
TrueTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from lFESpace to hFESpace...
virtual ~L2ProjectionH1Space()
virtual void Mult(const Vector &x, Vector &y) const override
Interpolation or prolongation of a true dof vector x to a true dof vector y.
virtual void SetRelTol(double p_rtol_)
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
virtual void SetAbsTol(double p_atol_)
virtual ~L2ProjectionGridTransfer()
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
L2ProjectionL2Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
virtual void SetRelTol(double p_rtol_)
virtual bool SupportsBackwardsOperator() const
TransferOperator(const FiniteElementSpace &lFESpace, const FiniteElementSpace &hFESpace)
Constructs a transfer operator from lFESpace to hFESpace.
Matrix-free transfer operator between finite element spaces working on true degrees of freedom...
virtual void MultTranspose(const Vector &x, Vector &y) const