xref: /petsc/src/mat/graphops/coarsen/impls/mis/mis.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
18be712e4SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
28be712e4SBarry Smith #include <../src/mat/impls/aij/seq/aij.h>
38be712e4SBarry Smith #include <../src/mat/impls/aij/mpi/mpiaij.h>
48be712e4SBarry Smith #include <petscsf.h>
58be712e4SBarry Smith 
68be712e4SBarry Smith #define MIS_NOT_DONE       -2
78be712e4SBarry Smith #define MIS_DELETED        -1
88be712e4SBarry Smith #define MIS_REMOVED        -3
98be712e4SBarry Smith #define MIS_IS_SELECTED(s) (s != MIS_DELETED && s != MIS_NOT_DONE && s != MIS_REMOVED)
108be712e4SBarry Smith 
118be712e4SBarry Smith /*
128be712e4SBarry Smith    MatCoarsenApply_MIS_private - parallel maximal independent set (MIS) with data locality info. MatAIJ specific!!!
138be712e4SBarry Smith 
148be712e4SBarry Smith    Input Parameter:
158be712e4SBarry Smith    . perm - serial permutation of rows of local to process in MIS
168be712e4SBarry Smith    . Gmat - global matrix of graph (data not defined)
178be712e4SBarry Smith    . strict_aggs - flag for whether to keep strict (non overlapping) aggregates in 'llist';
188be712e4SBarry Smith 
198be712e4SBarry Smith    Output Parameter:
208be712e4SBarry Smith    . a_selected - IS of selected vertices, includes 'ghost' nodes at end with natural local indices
218be712e4SBarry Smith    . a_locals_llist - array of list of nodes rooted at selected nodes
228be712e4SBarry Smith */
MatCoarsenApply_MIS_private(IS perm,Mat Gmat,PetscBool strict_aggs,PetscCoarsenData ** a_locals_llist)238be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MIS_private(IS perm, Mat Gmat, PetscBool strict_aggs, PetscCoarsenData **a_locals_llist)
248be712e4SBarry Smith {
258be712e4SBarry Smith   Mat_SeqAIJ       *matA, *matB = NULL;
268be712e4SBarry Smith   Mat_MPIAIJ       *mpimat = NULL;
278be712e4SBarry Smith   MPI_Comm          comm;
288be712e4SBarry Smith   PetscInt          num_fine_ghosts, kk, n, ix, j, *idx, *ii, Iend, my0, nremoved, gid, lid, cpid, lidj, sgid, t1, t2, slid, nDone, nselected = 0, state, statej;
298be712e4SBarry Smith   PetscInt         *cpcol_gid, *cpcol_state, *lid_cprowID, *lid_gid, *cpcol_sel_gid, *icpcol_gid, *lid_state, *lid_parent_gid = NULL, nrm_tot = 0;
308be712e4SBarry Smith   PetscBool        *lid_removed;
318be712e4SBarry Smith   PetscBool         isMPI, isAIJ, isOK;
328be712e4SBarry Smith   const PetscInt   *perm_ix;
338be712e4SBarry Smith   const PetscInt    nloc = Gmat->rmap->n;
348be712e4SBarry Smith   PetscCoarsenData *agg_lists;
358be712e4SBarry Smith   PetscLayout       layout;
368be712e4SBarry Smith   PetscSF           sf;
378be712e4SBarry Smith   IS                info_is;
388be712e4SBarry Smith 
398be712e4SBarry Smith   PetscFunctionBegin;
408be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
418be712e4SBarry Smith   PetscCall(ISCreate(comm, &info_is));
428be712e4SBarry Smith   PetscCall(PetscInfo(info_is, "mis: nloc = %d\n", (int)nloc));
438be712e4SBarry Smith   /* get submatrices */
448be712e4SBarry Smith   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPI));
458be712e4SBarry Smith   if (isMPI) {
468be712e4SBarry Smith     mpimat = (Mat_MPIAIJ *)Gmat->data;
478be712e4SBarry Smith     matA   = (Mat_SeqAIJ *)mpimat->A->data;
488be712e4SBarry Smith     matB   = (Mat_SeqAIJ *)mpimat->B->data;
498be712e4SBarry Smith     /* force compressed storage of B */
508be712e4SBarry Smith     PetscCall(MatCheckCompressedRow(mpimat->B, matB->nonzerorowcnt, &matB->compressedrow, matB->i, Gmat->rmap->n, -1.0));
518be712e4SBarry Smith   } else {
52c376f201SBarry Smith     matA = (Mat_SeqAIJ *)Gmat->data;
538be712e4SBarry Smith     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isAIJ));
548be712e4SBarry Smith     PetscCheck(isAIJ, comm, PETSC_ERR_PLIB, "Require AIJ matrix.");
558be712e4SBarry Smith   }
568be712e4SBarry Smith   PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend));
57c376f201SBarry Smith   PetscCall(PetscMalloc4(nloc, &lid_gid, nloc, &lid_cprowID, nloc, &lid_removed, nloc, &lid_state));
58c376f201SBarry Smith   if (strict_aggs) PetscCall(PetscMalloc1(nloc, &lid_parent_gid));
59c376f201SBarry Smith   if (isMPI) {
608be712e4SBarry Smith     for (kk = 0, gid = my0; kk < nloc; kk++, gid++) lid_gid[kk] = gid;
618be712e4SBarry Smith     PetscCall(VecGetLocalSize(mpimat->lvec, &num_fine_ghosts));
62c376f201SBarry Smith     PetscCall(PetscMalloc2(num_fine_ghosts, &cpcol_gid, num_fine_ghosts, &cpcol_state));
638be712e4SBarry Smith     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)Gmat), &sf));
648be712e4SBarry Smith     PetscCall(MatGetLayouts(Gmat, &layout, NULL));
658be712e4SBarry Smith     PetscCall(PetscSFSetGraphLayout(sf, layout, num_fine_ghosts, NULL, PETSC_COPY_VALUES, mpimat->garray));
668be712e4SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
678be712e4SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
688be712e4SBarry Smith     for (kk = 0; kk < num_fine_ghosts; kk++) cpcol_state[kk] = MIS_NOT_DONE;
698be712e4SBarry Smith   } else num_fine_ghosts = 0;
708be712e4SBarry Smith 
718be712e4SBarry Smith   /* has ghost nodes for !strict and uses local indexing (yuck) */
728be712e4SBarry Smith   PetscCall(PetscCDCreate(strict_aggs ? nloc : num_fine_ghosts + nloc, &agg_lists));
738be712e4SBarry Smith   if (a_locals_llist) *a_locals_llist = agg_lists;
748be712e4SBarry Smith 
758be712e4SBarry Smith   /* need an inverse map - locals */
768be712e4SBarry Smith   for (kk = 0; kk < nloc; kk++) {
778be712e4SBarry Smith     lid_cprowID[kk] = -1;
788be712e4SBarry Smith     lid_removed[kk] = PETSC_FALSE;
798be712e4SBarry Smith     if (strict_aggs) lid_parent_gid[kk] = -1.0;
808be712e4SBarry Smith     lid_state[kk] = MIS_NOT_DONE;
818be712e4SBarry Smith   }
828be712e4SBarry Smith   /* set index into cmpressed row 'lid_cprowID' */
838be712e4SBarry Smith   if (matB) {
848be712e4SBarry Smith     for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
858be712e4SBarry Smith       lid = matB->compressedrow.rindex[ix];
868be712e4SBarry Smith       if (lid >= 0) lid_cprowID[lid] = ix;
878be712e4SBarry Smith     }
888be712e4SBarry Smith   }
898be712e4SBarry Smith   /* MIS */
908be712e4SBarry Smith   nremoved = nDone = 0;
918be712e4SBarry Smith 
928be712e4SBarry Smith   PetscCall(ISGetIndices(perm, &perm_ix));
938be712e4SBarry Smith   while (nDone < nloc || PETSC_TRUE) { /* asynchronous not implemented */
948be712e4SBarry Smith     /* check all vertices */
958be712e4SBarry Smith     for (kk = 0; kk < nloc; kk++) {
968be712e4SBarry Smith       lid   = perm_ix[kk];
978be712e4SBarry Smith       state = lid_state[lid];
988be712e4SBarry Smith       if (lid_removed[lid]) continue;
998be712e4SBarry Smith       if (state == MIS_NOT_DONE) {
1008be712e4SBarry Smith         /* parallel test, delete if selected ghost */
1018be712e4SBarry Smith         isOK = PETSC_TRUE;
1028be712e4SBarry Smith         if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
1038be712e4SBarry Smith           ii  = matB->compressedrow.i;
1048be712e4SBarry Smith           n   = ii[ix + 1] - ii[ix];
1058be712e4SBarry Smith           idx = matB->j + ii[ix];
1068be712e4SBarry Smith           for (j = 0; j < n; j++) {
1078be712e4SBarry Smith             cpid   = idx[j]; /* compressed row ID in B mat */
1088be712e4SBarry Smith             gid    = cpcol_gid[cpid];
1098be712e4SBarry Smith             statej = cpcol_state[cpid];
110835f2295SStefano Zampini             PetscCheck(!MIS_IS_SELECTED(statej), PETSC_COMM_SELF, PETSC_ERR_SUP, "selected ghost: %" PetscInt_FMT, gid);
1118be712e4SBarry Smith             if (statej == MIS_NOT_DONE && gid >= Iend) { /* should be (pe>rank), use gid as pe proxy */
1128be712e4SBarry Smith               isOK = PETSC_FALSE;                        /* can not delete */
1138be712e4SBarry Smith               break;
1148be712e4SBarry Smith             }
1158be712e4SBarry Smith           }
1168be712e4SBarry Smith         } /* parallel test */
1178be712e4SBarry Smith         if (isOK) { /* select or remove this vertex */
1188be712e4SBarry Smith           nDone++;
1198be712e4SBarry Smith           /* check for singleton */
1208be712e4SBarry Smith           ii = matA->i;
1218be712e4SBarry Smith           n  = ii[lid + 1] - ii[lid];
1228be712e4SBarry Smith           if (n < 2) {
1238be712e4SBarry Smith             /* if I have any ghost adj then not a sing */
1248be712e4SBarry Smith             ix = lid_cprowID[lid];
1258be712e4SBarry Smith             if (ix == -1 || !(matB->compressedrow.i[ix + 1] - matB->compressedrow.i[ix])) {
1268be712e4SBarry Smith               nremoved++;
1278be712e4SBarry Smith               nrm_tot++;
1288be712e4SBarry Smith               lid_removed[lid] = PETSC_TRUE;
1298be712e4SBarry Smith               continue;
1308be712e4SBarry Smith               // lid_state[lidj] = MIS_REMOVED; /* add singleton to MIS (can cause low rank with elasticity on fine grid) */
1318be712e4SBarry Smith             }
1328be712e4SBarry Smith           }
1338be712e4SBarry Smith           /* SELECTED state encoded with global index */
1348be712e4SBarry Smith           lid_state[lid] = lid + my0;
1358be712e4SBarry Smith           nselected++;
1368be712e4SBarry Smith           if (strict_aggs) {
1378be712e4SBarry Smith             PetscCall(PetscCDAppendID(agg_lists, lid, lid + my0));
1388be712e4SBarry Smith           } else {
1398be712e4SBarry Smith             PetscCall(PetscCDAppendID(agg_lists, lid, lid));
1408be712e4SBarry Smith           }
1418be712e4SBarry Smith           /* delete local adj */
1428be712e4SBarry Smith           idx = matA->j + ii[lid];
1438be712e4SBarry Smith           for (j = 0; j < n; j++) {
1448be712e4SBarry Smith             lidj   = idx[j];
1458be712e4SBarry Smith             statej = lid_state[lidj];
1468be712e4SBarry Smith             if (statej == MIS_NOT_DONE) {
1478be712e4SBarry Smith               nDone++;
1488be712e4SBarry Smith               if (strict_aggs) {
1498be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lid, lidj + my0));
1508be712e4SBarry Smith               } else {
1518be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lid, lidj));
1528be712e4SBarry Smith               }
1538be712e4SBarry Smith               lid_state[lidj] = MIS_DELETED; /* delete this */
1548be712e4SBarry Smith             }
1558be712e4SBarry Smith           }
1568be712e4SBarry Smith           /* delete ghost adj of lid - deleted ghost done later for strict_aggs */
1578be712e4SBarry Smith           if (!strict_aggs) {
1588be712e4SBarry Smith             if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
1598be712e4SBarry Smith               ii  = matB->compressedrow.i;
1608be712e4SBarry Smith               n   = ii[ix + 1] - ii[ix];
1618be712e4SBarry Smith               idx = matB->j + ii[ix];
1628be712e4SBarry Smith               for (j = 0; j < n; j++) {
1638be712e4SBarry Smith                 cpid   = idx[j]; /* compressed row ID in B mat */
1648be712e4SBarry Smith                 statej = cpcol_state[cpid];
1658be712e4SBarry Smith                 if (statej == MIS_NOT_DONE) PetscCall(PetscCDAppendID(agg_lists, lid, nloc + cpid));
1668be712e4SBarry Smith               }
1678be712e4SBarry Smith             }
1688be712e4SBarry Smith           }
1698be712e4SBarry Smith         } /* selected */
1708be712e4SBarry Smith       } /* not done vertex */
1718be712e4SBarry Smith     } /* vertex loop */
1728be712e4SBarry Smith 
1738be712e4SBarry Smith     /* update ghost states and count todos */
174c376f201SBarry Smith     if (isMPI) {
1758be712e4SBarry Smith       /* scatter states, check for done */
1768be712e4SBarry Smith       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
1778be712e4SBarry Smith       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
1788be712e4SBarry Smith       ii = matB->compressedrow.i;
1798be712e4SBarry Smith       for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
1808be712e4SBarry Smith         lid   = matB->compressedrow.rindex[ix]; /* local boundary node */
1818be712e4SBarry Smith         state = lid_state[lid];
1828be712e4SBarry Smith         if (state == MIS_NOT_DONE) {
1838be712e4SBarry Smith           /* look at ghosts */
1848be712e4SBarry Smith           n   = ii[ix + 1] - ii[ix];
1858be712e4SBarry Smith           idx = matB->j + ii[ix];
1868be712e4SBarry Smith           for (j = 0; j < n; j++) {
1878be712e4SBarry Smith             cpid   = idx[j]; /* compressed row ID in B mat */
1888be712e4SBarry Smith             statej = cpcol_state[cpid];
1898be712e4SBarry Smith             if (MIS_IS_SELECTED(statej)) { /* lid is now deleted, do it */
1908be712e4SBarry Smith               nDone++;
1918be712e4SBarry Smith               lid_state[lid] = MIS_DELETED; /* delete this */
1928be712e4SBarry Smith               if (!strict_aggs) {
1938be712e4SBarry Smith                 lidj = nloc + cpid;
1948be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lidj, lid));
1958be712e4SBarry Smith               } else {
1968be712e4SBarry Smith                 sgid                = cpcol_gid[cpid];
1978be712e4SBarry Smith                 lid_parent_gid[lid] = sgid; /* keep track of proc that I belong to */
1988be712e4SBarry Smith               }
1998be712e4SBarry Smith               break;
2008be712e4SBarry Smith             }
2018be712e4SBarry Smith           }
2028be712e4SBarry Smith         }
2038be712e4SBarry Smith       }
2048be712e4SBarry Smith       /* all done? */
2058be712e4SBarry Smith       t1 = nloc - nDone;
206462c564dSBarry Smith       PetscCallMPI(MPIU_Allreduce(&t1, &t2, 1, MPIU_INT, MPI_SUM, comm)); /* synchronous version */
2078be712e4SBarry Smith       if (!t2) break;
2088be712e4SBarry Smith     } else break; /* all done */
2098be712e4SBarry Smith   } /* outer parallel MIS loop */
2108be712e4SBarry Smith   PetscCall(ISRestoreIndices(perm, &perm_ix));
2118be712e4SBarry Smith   PetscCall(PetscInfo(info_is, "\t removed %" PetscInt_FMT " of %" PetscInt_FMT " vertices.  %" PetscInt_FMT " selected.\n", nremoved, nloc, nselected));
2128be712e4SBarry Smith 
2138be712e4SBarry Smith   /* tell adj who my lid_parent_gid vertices belong to - fill in agg_lists selected ghost lists */
2148be712e4SBarry Smith   if (strict_aggs && matB) {
2158be712e4SBarry Smith     /* need to copy this to free buffer -- should do this globally */
216c376f201SBarry Smith     PetscCall(PetscMalloc2(num_fine_ghosts, &cpcol_sel_gid, num_fine_ghosts, &icpcol_gid));
2178be712e4SBarry Smith     for (cpid = 0; cpid < num_fine_ghosts; cpid++) icpcol_gid[cpid] = cpcol_gid[cpid];
2188be712e4SBarry Smith 
2198be712e4SBarry Smith     /* get proc of deleted ghost */
2208be712e4SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
2218be712e4SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
2228be712e4SBarry Smith     for (cpid = 0; cpid < num_fine_ghosts; cpid++) {
2238be712e4SBarry Smith       sgid = cpcol_sel_gid[cpid];
2248be712e4SBarry Smith       gid  = icpcol_gid[cpid];
2258be712e4SBarry Smith       if (sgid >= my0 && sgid < Iend) { /* I own this deleted */
2268be712e4SBarry Smith         slid = sgid - my0;
2278be712e4SBarry Smith         PetscCall(PetscCDAppendID(agg_lists, slid, gid));
2288be712e4SBarry Smith       }
2298be712e4SBarry Smith     }
230c376f201SBarry Smith     PetscCall(PetscFree2(cpcol_sel_gid, icpcol_gid));
2318be712e4SBarry Smith   }
232c376f201SBarry Smith   if (isMPI) {
2338be712e4SBarry Smith     PetscCall(PetscSFDestroy(&sf));
234c376f201SBarry Smith     PetscCall(PetscFree2(cpcol_gid, cpcol_state));
2358be712e4SBarry Smith   }
236c376f201SBarry Smith   PetscCall(PetscFree4(lid_gid, lid_cprowID, lid_removed, lid_state));
2378be712e4SBarry Smith   if (strict_aggs) {
2388be712e4SBarry Smith     // check sizes -- all vertices must get in graph
2398be712e4SBarry Smith     PetscInt aa[2] = {0, nrm_tot}, bb[2], MM;
240c376f201SBarry Smith 
241c376f201SBarry Smith     PetscCall(PetscFree(lid_parent_gid));
2428be712e4SBarry Smith     PetscCall(MatGetSize(Gmat, &MM, NULL));
2438be712e4SBarry Smith     // check sizes -- all vertices must get in graph
2448be712e4SBarry Smith     PetscCall(PetscCDCount(agg_lists, &aa[0]));
245462c564dSBarry Smith     PetscCallMPI(MPIU_Allreduce(aa, bb, 2, MPIU_INT, MPI_SUM, comm));
2468be712e4SBarry Smith     if (MM != bb[0]) PetscCall(PetscInfo(info_is, "Warning: N = %" PetscInt_FMT ", sum of aggregates %" PetscInt_FMT ", %" PetscInt_FMT " removed total\n", MM, bb[0], bb[1]));
247c376f201SBarry Smith     PetscCheck(MM >= bb[0], comm, PETSC_ERR_PLIB, "Sum of aggs is too large");
2488be712e4SBarry Smith   }
2498be712e4SBarry Smith   PetscCall(ISDestroy(&info_is));
2508be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2518be712e4SBarry Smith }
2528be712e4SBarry Smith 
2538be712e4SBarry Smith /*
2548be712e4SBarry Smith    MIS coarsen, simple greedy.
2558be712e4SBarry Smith */
MatCoarsenApply_MIS(MatCoarsen coarse)2568be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MIS(MatCoarsen coarse)
2578be712e4SBarry Smith {
2588be712e4SBarry Smith   Mat mat = coarse->graph;
2598be712e4SBarry Smith 
2608be712e4SBarry Smith   PetscFunctionBegin;
2618be712e4SBarry Smith   if (!coarse->perm) {
2628be712e4SBarry Smith     IS       perm;
2638be712e4SBarry Smith     PetscInt n, m;
2648be712e4SBarry Smith     MPI_Comm comm;
2658be712e4SBarry Smith 
2668be712e4SBarry Smith     PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
2678be712e4SBarry Smith     PetscCall(MatGetLocalSize(mat, &m, &n));
2688be712e4SBarry Smith     PetscCall(ISCreateStride(comm, m, 0, 1, &perm));
2698be712e4SBarry Smith     PetscCall(MatCoarsenApply_MIS_private(perm, mat, coarse->strict_aggs, &coarse->agg_lists));
2708be712e4SBarry Smith     PetscCall(ISDestroy(&perm));
2718be712e4SBarry Smith   } else {
2728be712e4SBarry Smith     PetscCall(MatCoarsenApply_MIS_private(coarse->perm, mat, coarse->strict_aggs, &coarse->agg_lists));
2738be712e4SBarry Smith   }
2748be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2758be712e4SBarry Smith }
2768be712e4SBarry Smith 
MatCoarsenView_MIS(MatCoarsen coarse,PetscViewer viewer)2778be712e4SBarry Smith static PetscErrorCode MatCoarsenView_MIS(MatCoarsen coarse, PetscViewer viewer)
2788be712e4SBarry Smith {
2798be712e4SBarry Smith   PetscMPIInt       rank;
280*9f196a02SMartin Diehl   PetscBool         isascii;
281e0b7e82fSBarry Smith   PetscViewerFormat format;
2828be712e4SBarry Smith 
2838be712e4SBarry Smith   PetscFunctionBegin;
2848be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank));
285*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
286e0b7e82fSBarry Smith   PetscCall(PetscViewerGetFormat(viewer, &format));
287*9f196a02SMartin Diehl   if (isascii && format == PETSC_VIEWER_ASCII_INFO_DETAIL) {
288bfa5bdfcSBarry Smith     if (coarse->agg_lists) {
2898be712e4SBarry Smith       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2908be712e4SBarry Smith       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] MIS aggregator\n", rank));
2918be712e4SBarry Smith       if (!rank) {
2928be712e4SBarry Smith         PetscCDIntNd *pos, *pos2;
2938be712e4SBarry Smith         for (PetscInt kk = 0; kk < coarse->agg_lists->size; kk++) {
2948be712e4SBarry Smith           PetscCall(PetscCDGetHeadPos(coarse->agg_lists, kk, &pos));
295835f2295SStefano Zampini           if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected %" PetscInt_FMT ": ", kk));
2968be712e4SBarry Smith           while (pos) {
2978be712e4SBarry Smith             PetscInt gid1;
2988be712e4SBarry Smith             PetscCall(PetscCDIntNdGetID(pos, &gid1));
2998be712e4SBarry Smith             PetscCall(PetscCDGetNextPos(coarse->agg_lists, kk, &pos));
300835f2295SStefano Zampini             PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", gid1));
3018be712e4SBarry Smith           }
3028be712e4SBarry Smith           if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
3038be712e4SBarry Smith         }
3048be712e4SBarry Smith       }
3058be712e4SBarry Smith       PetscCall(PetscViewerFlush(viewer));
3068be712e4SBarry Smith       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
307bfa5bdfcSBarry Smith     } else {
308bfa5bdfcSBarry Smith       PetscCall(PetscViewerASCIIPrintf(viewer, "  MIS aggregator lists are not available\n"));
309bfa5bdfcSBarry Smith     }
3108be712e4SBarry Smith   }
3118be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3128be712e4SBarry Smith }
3138be712e4SBarry Smith 
3148be712e4SBarry Smith /*MC
315c376f201SBarry Smith    MATCOARSENMIS - Creates a coarsening object that uses a maximal independent set (MIS) algorithm
3168be712e4SBarry Smith 
3178be712e4SBarry Smith    Collective
3188be712e4SBarry Smith 
3198be712e4SBarry Smith    Input Parameter:
3208be712e4SBarry Smith .  coarse - the coarsen context
3218be712e4SBarry Smith 
3228be712e4SBarry Smith    Level: beginner
3238be712e4SBarry Smith 
3248be712e4SBarry Smith .seealso: `MatCoarsen`, `MatCoarsenApply()`, `MatCoarsenGetData()`, `MatCoarsenSetType()`, `MatCoarsenType`
3258be712e4SBarry Smith M*/
MatCoarsenCreate_MIS(MatCoarsen coarse)3268be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_MIS(MatCoarsen coarse)
3278be712e4SBarry Smith {
3288be712e4SBarry Smith   PetscFunctionBegin;
3298be712e4SBarry Smith   coarse->ops->apply = MatCoarsenApply_MIS;
3308be712e4SBarry Smith   coarse->ops->view  = MatCoarsenView_MIS;
3318be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3328be712e4SBarry Smith }
333