14 #include "../../bilinearform.hpp"
15 #include "../../fespace.hpp"
16 #include "../../pfespace.hpp"
17 #include "../../../general/forall.hpp"
20 #include "../interface/restriction.hpp"
21 #include "../interface/ceed.hpp"
37 ConstrainedOperator(CeedOperator oper,
const Array<int> &ess_tdofs_,
40 ~ConstrainedOperator();
41 void Mult(
const Vector& x, Vector& y)
const;
42 CeedOperator GetCeedOperator()
const;
43 const Array<int> &GetEssentialTrueDofs()
const;
48 ceed::Operator *unconstrained_op;
54 const Array<int> &ess_tdofs_,
56 : ess_tdofs(ess_tdofs_), P(P_)
58 unconstrained_op =
new ceed::Operator(oper);
60 height = width = rap->
Height();
61 bool own_rap = (rap != unconstrained_op);
67 : ConstrainedOperator(oper, Array<int>(), P_)
72 delete constrained_op;
73 delete unconstrained_op;
78 constrained_op->Mult(x, y);
81 CeedOperator ConstrainedOperator::GetCeedOperator()
const
83 return unconstrained_op->GetCeedOperator();
86 const Array<int> &ConstrainedOperator::GetEssentialTrueDofs()
const
100 CeedSize in_len, out_len;
101 int ierr = CeedOperatorGetActiveVectorLengths(oper, &in_len, &out_len);
103 *size = (CeedInt)in_len;
104 MFEM_VERIFY(in_len == out_len,
"not a square CeedOperator");
105 MFEM_VERIFY(in_len == *size,
"size overflow");
112 CeedOperator ceed_op = op.GetCeedOperator();
113 const Array<int> &ess_tdofs = op.GetEssentialTrueDofs();
119 ierr = CeedVectorCreate(internal::ceed, length, &diagceed); PCeedChk(ierr);
121 ierr = CeedGetPreferredMemType(internal::ceed, &mem); PCeedChk(ierr);
126 Vector local_diag(length);
127 CeedScalar *ptr = (mem == CEED_MEM_HOST) ? local_diag.
HostWrite() :
128 local_diag.
Write(
true);
129 ierr = CeedVectorSetArray(diagceed, mem, CEED_USE_POINTER, ptr);
131 ierr = CeedOperatorLinearAssembleDiagonal(ceed_op, diagceed,
132 CEED_REQUEST_IMMEDIATE);
134 ierr = CeedVectorTakeArray(diagceed, mem, NULL); PCeedChk(ierr);
149 const int cheb_order = 3;
154 const double jacobi_scale = 0.65;
157 ierr = CeedVectorDestroy(&diagceed); PCeedChk(ierr);
164 class AssembledAMG :
public Solver
169 MFEM_ASSERT(P != NULL,
"Provided HypreParMatrix is invalid!");
170 height = width = oper.Height();
173 const Array<int> ess_tdofs = oper.GetEssentialTrueDofs();
181 op_assembled =
RAP(&hypre_local, P);
183 HypreParMatrix *mat_e = op_assembled->EliminateRowsCols(ess_tdofs);
185 amg =
new HypreBoomerAMG(*op_assembled);
186 amg->SetPrintLevel(0);
189 void Mult(
const Vector &x, Vector &y)
const override { amg->Mult(x, y); }
197 SparseMatrix *mat_local;
198 HypreParMatrix *op_assembled;
202 #endif // MFEM_USE_MPI
209 ho_boundary_ones = 0.0;
210 const int *ho_ess_tdofs_h = ho_ess_tdofs.
HostRead();
211 for (
int i=0; i<ho_ess_tdofs.
Size(); ++i)
213 ho_boundary_ones[ho_ess_tdofs_h[i]] = 1.0;
217 auto lobo = lo_boundary_ones.HostRead();
218 for (
int i = 0; i < lo_boundary_ones.Size(); ++i)
222 alg_lo_ess_tdofs.
Append(i);
235 MFEM_ABORT(
"This integrator does not support Ceed!");
243 ierr = CeedCompositeOperatorCreate(internal::ceed, &op); PCeedChk(ierr);
245 MFEM_VERIFY(form.
GetBBFI()->Size() == 0,
246 "Not implemented for this integrator!");
247 MFEM_VERIFY(form.
GetFBFI()->Size() == 0,
248 "Not implemented for this integrator!");
249 MFEM_VERIFY(form.
GetBFBFI()->Size() == 0,
250 "Not implemented for this integrator!");
254 for (
int i = 0; i < bffis->Size(); ++i)
263 CeedElemRestriction er,
270 ierr = CeedOperatorIsComposite(op, &isComposite); PCeedChk(ierr);
271 MFEM_ASSERT(isComposite,
"");
273 CeedOperator op_coarse;
274 ierr = CeedCompositeOperatorCreate(internal::ceed,
275 &op_coarse); PCeedChk(ierr);
278 CeedOperator *subops;
279 #if CEED_VERSION_GE(0, 10, 2)
280 ierr = CeedCompositeOperatorGetNumSub(op, &nsub); PCeedChk(ierr);
281 ierr = CeedCompositeOperatorGetSubList(op, &subops); PCeedChk(ierr);
283 ierr = CeedOperatorGetNumSub(op, &nsub); PCeedChk(ierr);
284 ierr = CeedOperatorGetSubList(op, &subops); PCeedChk(ierr);
286 for (
int isub=0; isub<nsub; ++isub)
288 CeedOperator subop = subops[isub];
289 CeedBasis basis_coarse, basis_c2f;
290 CeedOperator subop_coarse;
292 &basis_c2f, &subop_coarse); PCeedChk(ierr);
295 ierr = CeedBasisDestroy(&basis_coarse); PCeedChk(ierr);
296 ierr = CeedBasisDestroy(&basis_c2f); PCeedChk(ierr);
297 ierr = CeedCompositeOperatorAddSub(op_coarse, subop_coarse);
299 ierr = CeedOperatorDestroy(&subop_coarse); PCeedChk(ierr);
304 AlgebraicMultigrid::AlgebraicMultigrid(
311 ceed_operators.
SetSize(nlevels);
320 for (
int ilevel=nlevels-2; ilevel>=0; --ilevel)
333 for (
int ilevel=0; ilevel<nlevels; ++ilevel)
337 ConstrainedOperator *op =
new ConstrainedOperator(
357 if (P_mat) { smoother =
new AssembledAMG(*op, P_mat); }
373 int AlgebraicInterpolation::Initialize(
374 Ceed ceed, CeedBasis basisctof,
375 CeedElemRestriction erestrictu_coarse, CeedElemRestriction erestrictu_fine)
380 ierr = CeedElemRestrictionGetLVectorSize(erestrictu_coarse, &width);
382 ierr = CeedElemRestrictionGetLVectorSize(erestrictu_fine, &height);
386 const int bp3_ncompu = 1;
387 CeedQFunction l_qf_restrict, l_qf_prolong;
388 ierr = CeedQFunctionCreateIdentity(ceed, bp3_ncompu, CEED_EVAL_NONE,
389 CEED_EVAL_INTERP, &l_qf_restrict); CeedChk(ierr);
390 ierr = CeedQFunctionCreateIdentity(ceed, bp3_ncompu, CEED_EVAL_INTERP,
391 CEED_EVAL_NONE, &l_qf_prolong); CeedChk(ierr);
393 qf_restrict = l_qf_restrict;
394 qf_prolong = l_qf_prolong;
396 CeedVector c_fine_multiplicity;
397 ierr = CeedVectorCreate(ceed, height, &c_fine_multiplicity); CeedChk(ierr);
398 ierr = CeedVectorSetValue(c_fine_multiplicity, 0.0); CeedChk(ierr);
402 ierr = CeedOperatorCreate(ceed, qf_restrict, CEED_QFUNCTION_NONE,
403 CEED_QFUNCTION_NONE, &op_restrict); CeedChk(ierr);
404 ierr = CeedOperatorSetField(op_restrict,
"input", erestrictu_fine,
405 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); CeedChk(ierr);
406 ierr = CeedOperatorSetField(op_restrict,
"output", erestrictu_coarse,
407 basisctof, CEED_VECTOR_ACTIVE); CeedChk(ierr);
411 ierr = CeedOperatorCreate(ceed, qf_prolong, CEED_QFUNCTION_NONE,
412 CEED_QFUNCTION_NONE, &op_interp); CeedChk(ierr);
413 ierr = CeedOperatorSetField(op_interp,
"input", erestrictu_coarse,
414 basisctof, CEED_VECTOR_ACTIVE); CeedChk(ierr);
415 ierr = CeedOperatorSetField(op_interp,
"output", erestrictu_fine,
416 CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); CeedChk(ierr);
418 ierr = CeedElemRestrictionGetMultiplicity(erestrictu_fine,
419 c_fine_multiplicity); CeedChk(ierr);
420 ierr = CeedVectorCreate(ceed, height, &fine_multiplicity_r); CeedChk(ierr);
422 CeedScalar* fine_r_data;
423 const CeedScalar* fine_data;
424 ierr = CeedVectorGetArrayWrite(fine_multiplicity_r, CEED_MEM_HOST,
425 &fine_r_data); CeedChk(ierr);
426 ierr = CeedVectorGetArrayRead(c_fine_multiplicity, CEED_MEM_HOST,
427 &fine_data); CeedChk(ierr);
428 for (CeedSize i = 0; i <
height; ++i)
430 fine_r_data[i] = 1.0 / fine_data[i];
433 ierr = CeedVectorRestoreArray(fine_multiplicity_r, &fine_r_data); CeedChk(ierr);
434 ierr = CeedVectorRestoreArrayRead(c_fine_multiplicity, &fine_data);
436 ierr = CeedVectorDestroy(&c_fine_multiplicity); CeedChk(ierr);
438 ierr = CeedVectorCreate(ceed, height, &fine_work); CeedChk(ierr);
440 ierr = CeedVectorCreate(ceed, height, &v_); CeedChk(ierr);
441 ierr = CeedVectorCreate(ceed, width, &u_); CeedChk(ierr);
446 int AlgebraicInterpolation::Finalize()
450 ierr = CeedQFunctionDestroy(&qf_restrict); CeedChk(ierr);
451 ierr = CeedQFunctionDestroy(&qf_prolong); CeedChk(ierr);
452 ierr = CeedOperatorDestroy(&op_interp); CeedChk(ierr);
453 ierr = CeedOperatorDestroy(&op_restrict); CeedChk(ierr);
454 ierr = CeedVectorDestroy(&fine_multiplicity_r); CeedChk(ierr);
455 ierr = CeedVectorDestroy(&fine_work); CeedChk(ierr);
461 Ceed ceed, CeedBasis basisctof,
462 CeedElemRestriction erestrictu_coarse,
463 CeedElemRestriction erestrictu_fine)
466 CeedSize lo_nldofs, ho_nldofs;
467 ierr = CeedElemRestrictionGetLVectorSize(erestrictu_coarse, &lo_nldofs);
469 ierr = CeedElemRestrictionGetLVectorSize(erestrictu_fine,
470 &ho_nldofs); PCeedChk(ierr);
471 height = (int)ho_nldofs;
472 width = (int)lo_nldofs;
473 MFEM_VERIFY(ho_nldofs == height,
"height overflow");
474 MFEM_VERIFY(lo_nldofs == width,
"width overflow");
476 ierr = Initialize(ceed, basisctof, erestrictu_coarse, erestrictu_fine);
483 ierr = CeedVectorDestroy(&v_); PCeedChk(ierr);
484 ierr = CeedVectorDestroy(&u_); PCeedChk(ierr);
487 ierr = CeedBasisDestroy(&basisctof_); PCeedChk(ierr);
498 CeedVectorGetCeed(a, &ceed);
500 CeedSize length, length2;
501 ierr = CeedVectorGetLength(a, &length); CeedChk(ierr);
502 ierr = CeedVectorGetLength(b, &length2); CeedChk(ierr);
503 if (length != length2)
505 return CeedError(ceed, 1,
"Vector sizes don't match");
511 mem = CEED_MEM_DEVICE;
518 const CeedScalar *b_data;
519 ierr = CeedVectorGetArray(a, mem, &a_data); CeedChk(ierr);
520 ierr = CeedVectorGetArrayRead(b, mem, &b_data); CeedChk(ierr);
521 MFEM_VERIFY(
int(length) == length,
"length overflow");
522 MFEM_FORALL(i, length,
523 {a_data[i] *= b_data[i];});
525 ierr = CeedVectorRestoreArray(a, &a_data); CeedChk(ierr);
526 ierr = CeedVectorRestoreArrayRead(b, &b_data); CeedChk(ierr);
534 const CeedScalar *in_ptr;
537 ierr = CeedGetPreferredMemType(internal::ceed, &mem); PCeedChk(ierr);
549 ierr = CeedVectorSetArray(u_, mem, CEED_USE_POINTER,
550 const_cast<CeedScalar*>(in_ptr)); PCeedChk(ierr);
551 ierr = CeedVectorSetArray(v_, mem, CEED_USE_POINTER,
552 out_ptr); PCeedChk(ierr);
554 ierr = CeedOperatorApply(op_interp, u_, v_,
555 CEED_REQUEST_IMMEDIATE); PCeedChk(ierr);
558 ierr = CeedVectorTakeArray(u_, mem, const_cast<CeedScalar**>(&in_ptr));
560 ierr = CeedVectorTakeArray(v_, mem, &out_ptr); PCeedChk(ierr);
568 ierr = CeedGetPreferredMemType(internal::ceed, &mem); PCeedChk(ierr);
569 const CeedScalar *in_ptr;
582 ierr = CeedVectorSetArray(v_, mem, CEED_USE_POINTER,
583 const_cast<CeedScalar*>(in_ptr)); PCeedChk(ierr);
584 ierr = CeedVectorSetArray(u_, mem, CEED_USE_POINTER,
585 out_ptr); PCeedChk(ierr);
588 ierr = CeedVectorGetLength(v_, &length); PCeedChk(ierr);
590 const CeedScalar *multiplicitydata;
591 CeedScalar *workdata;
592 ierr = CeedVectorGetArrayRead(fine_multiplicity_r, mem,
593 &multiplicitydata); PCeedChk(ierr);
594 ierr = CeedVectorGetArrayWrite(fine_work, mem, &workdata); PCeedChk(ierr);
595 MFEM_VERIFY((
int)length == length,
"length overflow");
596 MFEM_FORALL(i, length,
597 {workdata[i] = in_ptr[i] * multiplicitydata[i];});
598 ierr = CeedVectorRestoreArrayRead(fine_multiplicity_r,
600 ierr = CeedVectorRestoreArray(fine_work, &workdata); PCeedChk(ierr);
602 ierr = CeedOperatorApply(op_restrict, fine_work, u_,
603 CEED_REQUEST_IMMEDIATE); PCeedChk(ierr);
605 ierr = CeedVectorTakeArray(v_, mem, const_cast<CeedScalar**>(&in_ptr));
607 ierr = CeedVectorTakeArray(u_, mem, &out_ptr); PCeedChk(ierr);
614 int current_order = order;
615 while (current_order > 0)
618 current_order = current_order/2;
633 ceed_interpolations.SetSize(nlevels-1);
634 R_tr.SetSize(nlevels-1);
638 current_order = order;
640 Ceed ceed = internal::ceed;
642 CeedElemRestriction er = fine_er;
654 for (
int ilevel=nlevels-2; ilevel>=0; --ilevel)
656 const int order_reduction = current_order - (current_order/2);
663 *
fespaces[ilevel+1], er, current_order, dim, order_reduction, gc);
671 *
fespaces[ilevel+1], er, current_order, dim, order_reduction);
673 current_order = current_order/2;
677 space->GetCeedCoarseToFine(),
678 space->GetCeedElemRestriction(),
690 prolongations[ilevel] = ceed_interpolations[ilevel]->SetupRAP(
691 space->GetProlongationMatrix(), R_tr[ilevel]);
695 er = space->GetCeedElemRestriction();
701 CeedElemRestriction fine_er,
705 ) : order_reduction(order_reduction_)
720 MFEM_VERIFY(
ndofs == ndofs_,
"ndofs overflow");
737 CeedElemRestriction fine_er,
740 int order_reduction_,
748 MFEM_VERIFY((
int)lsize == lsize,
"size overflow");
755 group_ldof.
MakeI(group_ldof_fine.
Size());
756 for (
int g=1; g<group_ldof_fine.
Size(); ++g)
758 int nldof_fine_g = group_ldof_fine.
RowSize(g);
759 const int *ldof_fine_g = group_ldof_fine.
GetRow(g);
760 for (
int i=0; i<nldof_fine_g; ++i)
762 int icoarse =
dof_map[ldof_fine_g[i]];
766 ldof_group[icoarse] = g;
771 for (
int g=1; g<group_ldof_fine.
Size(); ++g)
773 int nldof_fine_g = group_ldof_fine.
RowSize(g);
774 const int *ldof_fine_g = group_ldof_fine.
GetRow(g);
775 for (
int i=0; i<nldof_fine_g; ++i)
777 int icoarse =
dof_map[ldof_fine_g[i]];
789 for (
int i=0; i<lsize; ++i)
791 int g = ldof_group[i];
794 ldof_ltdof[i] = ltsize;
799 gc->
Bcast(ldof_ltdof);
802 for (
int j=0; j<lsize; ++j)
806 int i = ldof_ltdof[j];
825 if (P_mat) {
return P_mat; }
828 MFEM_VERIFY(pmesh != NULL,
"");
832 int ltsize = P->
Width();
836 MPI_Comm comm = pmesh->
GetComm();
840 if (HYPRE_AssumedPartitionCheck())
844 MPI_Request *requests =
new MPI_Request[2*nsize];
845 MPI_Status *statuses =
new MPI_Status[2*nsize];
846 tdof_nb_offsets.
SetSize(nsize+1);
847 tdof_nb_offsets[0] = tdof_offsets[0];
850 int request_counter = 0;
851 for (
int i = 1; i <= nsize; i++)
853 MPI_Irecv(&tdof_nb_offsets[i], 1, HYPRE_MPI_INT,
855 &requests[request_counter++]);
857 for (
int i = 1; i <= nsize; i++)
859 MPI_Isend(&tdof_nb_offsets[0], 1, HYPRE_MPI_INT,
861 &requests[request_counter++]);
863 MPI_Waitall(request_counter, requests, statuses);
884 i_diag[0] = i_offd[0] = 0;
885 diag_counter = offd_counter = 0;
886 for (
int i_ldof = 0; i_ldof < lsize; i_ldof++)
888 int g = ldof_group[i_ldof];
889 int i_ltdof = ldof_ltdof[i_ldof];
892 j_diag[diag_counter++] = i_ltdof;
897 if (HYPRE_AssumedPartitionCheck())
908 cmap_j_offd[offd_counter].one = global_tdof_number;
909 cmap_j_offd[offd_counter].two = offd_counter;
912 i_diag[i_ldof+1] = diag_counter;
913 i_offd[i_ldof+1] = offd_counter;
916 SortPairs<HYPRE_BigInt, int>(cmap_j_offd, offd_counter);
918 for (
int i = 0; i < offd_counter; i++)
920 cmap[i] = cmap_j_offd[i].one;
921 j_offd[cmap_j_offd[i].two] = i;
926 row_starts, col_starts,
927 i_diag, j_diag, i_offd, j_offd,
945 #endif // MFEM_USE_MPI
947 #endif // MFEM_USE_CEED
953 "AlgebraicSolver requires a Ceed device");
956 "AlgebraicSolver requires partial assembly or fully matrix-free.");
958 "AlgebraicSolver requires tensor product basis functions.");
963 MFEM_ABORT(
"AlgebraicSolver requires Ceed support");
978 multigrid->
Mult(x, y);
int GetGroupMasterRank(int g) const
Return the rank of the group master for group 'g'.
int Size() const
Return the logical size of the array.
int CeedATPMGOperator(CeedOperator oper, int order_reduction, CeedElemRestriction coarse_er, CeedBasis coarse_basis_in, CeedBasis basis_ctof_in, CeedOperator *out)
ConstrainedOperator(Operator *A, const Array< int > &list, bool own_A=false, DiagonalPolicy diag_policy=DIAG_ONE)
Constructor from a general Operator and a list of essential indices/dofs.
MPI_Comm GetComm() const
MPI communicator.
int ndofs
Number of degrees of freedom. Number of unknowns is ndofs * vdim.
void Bcast(T *ldata, int layout) const
Broadcast within each group where the master is the root.
Chebyshev accelerated smoothing with given vector, no matrix necessary.
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
CeedElemRestriction GetCeedElemRestriction() const
AlgebraicInterpolation(Ceed ceed, CeedBasis basisctof, CeedElemRestriction erestrictu_coarse, CeedElemRestriction erestrictu_fine)
virtual void Finalize(int skip_zeros=1)
Finalize the matrix initialization, switching the storage format from LIL to CSR. ...
void SetSize(int s)
Resize the vector to size s.
Geometric multigrid associated with a hierarchy of finite element spaces.
const T * HostRead() const
Shortcut for mfem::Read(a.GetMemory(), a.Size(), false).
void Mult(const Table &A, const Table &B, Table &C)
C = A * B (as boolean matrices)
int Width() const
Get the width (size of input) of the Operator. Synonym with NumCols().
virtual void Mult(const Vector &x, Vector &y) const
Constrained operator action.
virtual double * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Abstract parallel finite element space.
bool UsesTensorBasis(const FiniteElementSpace &fes)
Return true if the mesh contains only one topology and the elements are tensor elements.
Extension of Multigrid object to algebraically generated coarse spaces.
GroupCommunicator * GetGroupCommunicator() const
Memory< double > & GetMemory()
Return a reference to the Memory object used by the Vector.
A way to use algebraic levels in a Multigrid object.
Array< Operator * > prolongations
Solver * BuildSmootherFromCeed(ConstrainedOperator &op, bool chebyshev)
virtual 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 ...
bool IAmMaster(int g) const
Return true if I am master for group 'g'.
CeedOperator & GetCeedOperator()
int CeedVectorPointwiseMult(CeedVector a, const CeedVector b)
Array< Array< int > * > essentialTrueDofs
virtual const Operator * GetProlongation() const
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator...
void Set(const int i, const int j, const double val)
Communicator performing operations within groups defined by a GroupTopology with arbitrary-size data ...
~AlgebraicInterpolation()
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Array< bool > ownedProlongations
int Append(const T &el)
Append element 'el' to array, resize if necessary.
virtual ~ConstrainedOperator()
Destructor: destroys the unconstrained Operator, if owned.
Mesh * GetMesh() const
Returns the mesh.
AlgebraicSolver(BilinearForm &form, const Array< int > &ess_tdofs)
Constructs algebraic multigrid hierarchy and solver.
virtual double * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
int GetNumNeighbors() const
Return the number of neighbors including the local processor.
Jacobi smoothing for a given bilinear form (no matrix necessary).
void Finalize()
Allocate internal buffers after the GroupLDofTable is defined.
void AddConnection(int r, int c)
void Mult(const Vector &x, Vector &y) const
Operator application: y=A(x).
HypreParMatrix * Dof_TrueDof_Matrix() const
The true dof-to-dof interpolation matrix.
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 ...
virtual const Operator * GetProlongationMatrix() const
The returned Operator is owned by the FiniteElementSpace.
Array< FiniteElementSpace * > fespaces
int GetNumLevels() const
Returns the number of levels in the hierarchy.
CeedBasis GetCeedCoarseToFine() const
int CeedATPMGElemRestriction(int order, int order_reduction, CeedElemRestriction er_in, CeedElemRestriction *er_out, CeedInt *&dof_map)
Take given (high-order) CeedElemRestriction and make a new CeedElemRestriction, which corresponds to ...
int Size() const
Returns the number of TYPE I elements.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
Multigrid interpolation operator in Ceed framework.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
const FiniteElementSpaceHierarchy & fespaces
int GetGroupMaster(int g) const
Return the neighbor index of the group master for a given group. Neighbor 0 is the local processor...
virtual const FiniteElementSpace & GetFESpaceAtLevel(int level) const
Returns the finite element space at the given level.
const GroupTopology & GetGroupTopology()
Get a reference to the associated GroupTopology object.
virtual void Mult(const mfem::Vector &x, mfem::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.
int GetNeighborRank(int i) const
Return the MPI rank of neighbor 'i'.
Array< bool > ownedMeshes
void AddAColumnInRow(int r)
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
GroupCommunicator & GroupComm()
Return a reference to the internal GroupCommunicator (on VDofs)
CeedOperator CoarsenCeedCompositeOperator(CeedOperator op, CeedElemRestriction er, CeedBasis c2f, int order_reduction)
bool DeviceCanUseCeed()
Function that determines if a CEED kernel should be used, based on the current mfem::Device configura...
int GetOrder(int i) const
Returns the polynomial degree of the i'th finite element.
~ParAlgebraicCoarseSpace()
CeedOperator CreateCeedCompositeOperatorFromBilinearForm(BilinearForm &form)
int CeedBasisATPMGCoarseToFine(Ceed ceed, int P1d, int dim, int order_reduction, CeedBasis *basisc2f)
Create coarse-to-fine basis, given number of input nodes and order reduction.
int CeedOperatorGetSize(CeedOperator oper, CeedInt *size)
int GetOrderReduction() const
void SetLTDofTable(const Array< int > &ldof_ltdof)
Initialize the internal group_ltdof Table.
static bool Allows(unsigned long b_mask)
Return true if any of the backends in the backend mask, b_mask, are allowed.
The transpose of a given operator. Switches the roles of the methods Mult() and MultTranspose().
int CeedOperatorFullAssemble(CeedOperator op, SparseMatrix **mat)
Assembles a CeedOperator as an mfem::SparseMatrix.
Table & GroupLDofTable()
Fill-in the returned Table reference to initialize the GroupCommunicator then call Finalize()...
AlgebraicCoarseSpace(FiniteElementSpace &fine_fes, CeedElemRestriction fine_er, int order, int dim, int order_reduction_)
int height
Dimension of the output / number of rows in the matrix.
void AddLevel(Operator *opr, Solver *smoother, bool ownOperator, bool ownSmoother)
Adds a level to the multigrid operator hierarchy.
CeedElemRestriction ceed_elem_restriction
Parallel version of AlgebraicCoarseSpace.
virtual void Mult(const Vector &x, Vector &y) const override
Application of the multigrid as a preconditioner.
virtual double * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
Mesh * mesh
The mesh that FE space lives on (not owned).
HYPRE_BigInt GetGlobalNumRows() const
Return the global number of rows.
HYPRE_BigInt * RowPart()
Returns the row partitioning.
void NewMemoryAndSize(const Memory< double > &mem, int s, bool own_mem)
Reset the Vector to use the given external Memory mem and size s.
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.
Biwise-OR of all device backends.
void AddToCompositeOperator(BilinearFormIntegrator *integ, CeedOperator op)
Square Operator for imposing essential boundary conditions using only the action, Mult()...
Operator * GetProlongationAtLevel(int level) const
Returns the prolongation operator from the finite element space at level to the finite element space ...
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Wrapper for hypre's ParCSR matrix class.
Class for parallel meshes.
[device] CUDA backend. Enabled when MFEM_USE_CUDA = YES.
virtual const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
virtual double * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
void GenerateOffsets(int N, HYPRE_BigInt loc_sizes[], Array< HYPRE_BigInt > *offsets[]) const
int width
Dimension of the input / number of columns in the matrix.
void InitRestriction(const FiniteElementSpace &fes, Ceed ceed, CeedElemRestriction *restr)
Initialize a CeedElemRestriction for non-mixed meshes.
void CoarsenEssentialDofs(const mfem::Operator &interp, const Array< int > &ho_ess_tdofs, Array< int > &alg_lo_ess_tdofs)
HypreParMatrix * GetProlongationHypreParMatrix()
Hierarchy of AlgebraicCoarseSpace objects for use in Multigrid object.
AlgebraicCoarseSpace & GetAlgebraicCoarseSpace(int level)