12 #ifndef MFEM_CEED_ALGEBRAIC_HPP 13 #define MFEM_CEED_ALGEBRAIC_HPP 15 #include "../../fespacehierarchy.hpp" 16 #include "../../multigrid.hpp" 17 #include "../interface/operator.hpp" 18 #include "../interface/ceed.hpp" 35 int order,
int dim,
int order_reduction_);
61 CeedElemRestriction fine_er,
81 #endif // MFEM_USE_MPI 92 Ceed ceed, CeedBasis basisctof,
93 CeedElemRestriction erestrictu_coarse,
94 CeedElemRestriction erestrictu_fine);
104 int Initialize(Ceed ceed, CeedBasis basisctof,
105 CeedElemRestriction erestrictu_coarse,
106 CeedElemRestriction erestrictu_fine);
109 CeedBasis basisctof_;
114 CeedQFunction qf_restrict, qf_prolong;
115 CeedOperator op_interp, op_restrict;
116 CeedVector fine_multiplicity_r;
117 CeedVector fine_work;
137 for (
int i=0; i<R_tr.Size(); ++i)
141 for (
int i=0; i<ceed_interpolations.Size(); ++i)
143 delete ceed_interpolations[i];
148 CeedElemRestriction fine_er;
178 #endif // MFEM_USE_CEED 211 #endif // MFEM_CEED_ALGEBRAIC_HPP AlgebraicInterpolation(Ceed ceed, CeedBasis basisctof, CeedElemRestriction erestrictu_coarse, CeedElemRestriction erestrictu_fine)
virtual const SparseMatrix * GetRestrictionMatrix() const override
The returned SparseMatrix is owned by the FiniteElementSpace.
Geometric multigrid associated with a hierarchy of finite element spaces.
Pointer to an Operator of a specified type.
virtual const mfem::Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
~AlgebraicSpaceHierarchy()
Extension of Multigrid object to algebraically generated coarse spaces.
A way to use algebraic levels in a Multigrid object.
GroupCommunicator * GetGroupCommunicator() const
CeedElemRestriction GetCeedElemRestriction() const
int GetNumLevels() const
Returns the number of levels in the hierarchy.
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
~AlgebraicInterpolation()
AlgebraicSolver(BilinearForm &form, const Array< int > &ess_tdofs)
Constructs algebraic multigrid hierarchy and solver.
virtual const mfem::Operator * GetProlongationMatrix() const override
The returned Operator is owned by the FiniteElementSpace.
Array< FiniteElementSpace * > fespaces
Multigrid interpolation operator in Ceed framework.
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
void SetOperator(const mfem::Operator &op)
Set/update the solver for the given operator.
virtual void MultTranspose(const mfem::Vector &x, mfem::Vector &y) const
Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an ...
~ParAlgebraicCoarseSpace()
virtual const SparseMatrix * GetRestrictionMatrix() const override
The returned SparseMatrix is owned by the FiniteElementSpace.
AlgebraicCoarseSpace(FiniteElementSpace &fine_fes, CeedElemRestriction fine_er, int order, int dim, int order_reduction_)
CeedElemRestriction ceed_elem_restriction
Wrapper for AlgebraicMultigrid object.
Parallel version of AlgebraicCoarseSpace.
int GetOrderReduction() const
CeedBasis GetCeedCoarseToFine() const
virtual void Mult(const mfem::Vector &x, mfem::Vector &y) const
Operator application: y=A(x).
virtual void SetOperator(const mfem::Operator &op) override
Not supported for multigrid.
Operator * SetupRAP(const Operator *Pi, const Operator *Po)
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P"...
ParAlgebraicCoarseSpace(FiniteElementSpace &fine_fes, CeedElemRestriction fine_er, int order, int dim, int order_reduction_, GroupCommunicator *gc_fine)
AlgebraicSpaceHierarchy(FiniteElementSpace &fespace)
Construct hierarchy based on finest FiniteElementSpace.
Wrapper for hypre's ParCSR matrix class.
AlgebraicMultigrid(AlgebraicSpaceHierarchy &hierarchy, BilinearForm &form, const Array< int > &ess_tdofs)
Constructs multigrid solver based on existing space hierarchy.
HypreParMatrix * GetProlongationHypreParMatrix()
Hierarchy of AlgebraicCoarseSpace objects for use in Multigrid object.
AlgebraicCoarseSpace & GetAlgebraicCoarseSpace(int level)