12#ifndef MFEM_TRANSFER_HPP
13#define MFEM_TRANSFER_HPP
144 bool own_mass_integ_ =
true);
311 std::pair<std::unique_ptr<SparseMatrix>,
341 std::unique_ptr<SparseMatrix>
AllocR();
348 std::unique_ptr<Operator>
R;
350 std::unique_ptr<Operator>
M_LH;
381 bool force_l2_space_ =
false)
470 const Operator* elem_restrict_lex_l;
471 const Operator* elem_restrict_lex_h;
Conjugate gradient method.
Data type dense matrix using column-major storage.
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
Abstract class for all finite elements.
Base class for transfer algorithms that construct transfer Operators between two finite element (FE) ...
virtual ~GridTransfer()
Virtual destructor.
FiniteElementSpace & dom_fes
Domain FE space.
OperatorHandle fw_t_oper
Forward true-dof operator.
FiniteElementSpace & ran_fes
Range FE space.
virtual const Operator & TrueBackwardOperator()
Return an Operator that transfers true-dof Vectors from the range FE space back to true-dof Vectors i...
OperatorHandle bw_t_oper
Backward true-dof operator.
virtual bool SupportsBackwardsOperator() const
virtual const Operator & ForwardOperator()=0
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
void SetOperatorType(Operator::Type type)
Set the desired Operator::Type for the construction of all operators defined by the underlying transf...
virtual const Operator & BackwardOperator()=0
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Operator::Type oper_type
Desired Operator::Type for the construction of all operators defined by the underlying transfer algor...
virtual const Operator & TrueForwardOperator()
Return an Operator that transfers true-dof Vectors from the domain FE space to true-dof Vectors in th...
const Operator & MakeTrueOperator(FiniteElementSpace &fes_in, FiniteElementSpace &fes_out, const Operator &oper, OperatorHandle &t_oper)
GridTransfer(FiniteElementSpace &dom_fes_, FiniteElementSpace &ran_fes_)
Transfer data between a coarse mesh and an embedded refined mesh using interpolation.
virtual ~InterpolationGridTransfer()
InterpolationGridTransfer(FiniteElementSpace &coarse_fes, FiniteElementSpace &fine_fes)
BilinearFormIntegrator * mass_integ
Ownership depends on own_mass_integ.
OperatorHandle F
Forward, coarse-to-fine, operator.
bool own_mass_integ
Ownership flag for mass_integ.
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
OperatorHandle B
Backward, fine-to-coarse, 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,...
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
void SetFromTDofsTranspose(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dual field coefficients given a vector of dual field coefficients on the tdofs and a finite elem...
virtual void Prolongate(const Vector &x, Vector &y) const
void TDofsListByVDim(const FiniteElementSpace &fes, int vdim, Array< int > &vdofs_list) const
Fills the vdofs_list array with a list of vdofs for a given vdim and a given finite element space.
std::unique_ptr< Solver > precon
void SetupPCG()
Sets up the PCG solver (sets parameters, operator, and preconditioner)
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
void GetTDofs(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers vector of tdofs given a vector of dofs and a finite element space.
void SetFromTDofs(const FiniteElementSpace &fes, const Vector &X, Vector &x) const
Sets dof values given a vector of tdofs and a finite element space.
virtual void SetAbsTol(real_t p_atol_)
Sets absolute tolerance in preconditioned conjugate gradient solver.
std::unique_ptr< Operator > RTxM_LH
L2ProjectionH1Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
void GetTDofsTranspose(const FiniteElementSpace &fes, const Vector &x, Vector &X) const
Recovers a vector of dual field coefficients on the tdofs given a vector of dual coefficients and a f...
virtual void Mult(const Vector &x, Vector &y) const
std::unique_ptr< SparseMatrix > AllocR()
Computes sparsity pattern and initializes R matrix.
std::unique_ptr< Operator > M_LH
void LumpedMassInverse(Vector &ML_inv) const
Returns the inverse of an on-rank lumped mass matrix.
std::pair< std::unique_ptr< SparseMatrix >, std::unique_ptr< SparseMatrix > > ComputeSparseRAndM_LH()
Computes on-rank R and M_LH matrices.
std::unique_ptr< Operator > R
virtual void MultTranspose(const Vector &x, Vector &y) const
virtual void SetRelTol(real_t p_rtol_)
Sets relative tolerance in preconditioned conjugate gradient solver.
virtual void SetRelTol(real_t p_rtol_)
No-op.
L2ProjectionL2Space(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
virtual void Prolongate(const Vector &x, Vector &y) const
virtual void MultTranspose(const Vector &x, Vector &y) const
virtual void Mult(const Vector &x, Vector &y) const
virtual void SetAbsTol(real_t p_atol_)
No-op.
virtual void ProlongateTranspose(const Vector &x, Vector &y) const
void BuildHo2Lor(int nel_ho, int nel_lor, const CoarseFineTransformations &cf_tr)
const FiniteElementSpace & fes_ho
virtual void ProlongateTranspose(const Vector &x, Vector &y) const =0
virtual void SetRelTol(real_t p_rtol_)=0
Sets relative tolerance in preconditioned conjugate gradient solver.
const FiniteElementSpace & fes_lor
void ElemMixedMass(Geometry::Type geom, const FiniteElement &fe_ho, const FiniteElement &fe_lor, ElementTransformation *el_tr, IntegrationPointTransformation &ip_tr, DenseMatrix &M_mixed_el) const
L2Projection(const FiniteElementSpace &fes_ho_, const FiniteElementSpace &fes_lor_)
virtual void SetAbsTol(real_t p_atol_)=0
Sets absolute tolerance in preconditioned conjugate gradient solver.
virtual void Prolongate(const Vector &x, Vector &y) const =0
L2Prolongation(const L2Projection &l2proj_)
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 ~L2Prolongation()
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Transfer data in L2 and H1 finite element spaces between a coarse mesh and an embedded refined mesh u...
L2Projection * F
Forward, coarse-to-fine, operator.
virtual bool SupportsBackwardsOperator() const
L2ProjectionGridTransfer(FiniteElementSpace &coarse_fes_, FiniteElementSpace &fine_fes_, bool force_l2_space_=false)
virtual const Operator & ForwardOperator()
Return an Operator that transfers GridFunctions from the domain FE space to GridFunctions in the rang...
L2Prolongation * B
Backward, fine-to-coarse, operator.
virtual ~L2ProjectionGridTransfer()
virtual const Operator & BackwardOperator()
Return an Operator that transfers GridFunctions from the range FE space back to GridFunctions in the ...
Pointer to an Operator of a specified type.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
Type
Enumeration defining IDs for some classes derived from Operator.
Matrix-free transfer operator between finite element spaces on the same mesh.
PRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~PRefinementTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
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...
Abstract parallel finite element space.
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...
TensorProductPRefinementTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator from lFESpace to hFESpace which have different FE collections.
virtual ~TensorProductPRefinementTransferOperator()
Destructor.
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.
virtual ~TransferOperator()
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 void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
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.
~TrueTransferOperator()
Destructor.
virtual void MultTranspose(const Vector &x, Vector &y) const override
Restriction by applying the transpose of the Mult method.
TrueTransferOperator(const FiniteElementSpace &lFESpace_, const FiniteElementSpace &hFESpace_)
Constructs a transfer operator working on true degrees of freedom from lFESpace to hFESpace.
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.