xref: /petsc/src/mat/graphops/coarsen/impls/mis/mis.c (revision 8be712e46db5d855f641c6bd97b4543e0efe65bd)
1*8be712e4SBarry Smith #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2*8be712e4SBarry Smith #include <../src/mat/impls/aij/seq/aij.h>
3*8be712e4SBarry Smith #include <../src/mat/impls/aij/mpi/mpiaij.h>
4*8be712e4SBarry Smith #include <petscsf.h>
5*8be712e4SBarry Smith 
6*8be712e4SBarry Smith #define MIS_NOT_DONE       -2
7*8be712e4SBarry Smith #define MIS_DELETED        -1
8*8be712e4SBarry Smith #define MIS_REMOVED        -3
9*8be712e4SBarry Smith #define MIS_IS_SELECTED(s) (s != MIS_DELETED && s != MIS_NOT_DONE && s != MIS_REMOVED)
10*8be712e4SBarry Smith 
11*8be712e4SBarry Smith /*
12*8be712e4SBarry Smith    MatCoarsenApply_MIS_private - parallel maximal independent set (MIS) with data locality info. MatAIJ specific!!!
13*8be712e4SBarry Smith 
14*8be712e4SBarry Smith    Input Parameter:
15*8be712e4SBarry Smith    . perm - serial permutation of rows of local to process in MIS
16*8be712e4SBarry Smith    . Gmat - global matrix of graph (data not defined)
17*8be712e4SBarry Smith    . strict_aggs - flag for whether to keep strict (non overlapping) aggregates in 'llist';
18*8be712e4SBarry Smith 
19*8be712e4SBarry Smith    Output Parameter:
20*8be712e4SBarry Smith    . a_selected - IS of selected vertices, includes 'ghost' nodes at end with natural local indices
21*8be712e4SBarry Smith    . a_locals_llist - array of list of nodes rooted at selected nodes
22*8be712e4SBarry Smith */
23*8be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MIS_private(IS perm, Mat Gmat, PetscBool strict_aggs, PetscCoarsenData **a_locals_llist)
24*8be712e4SBarry Smith {
25*8be712e4SBarry Smith   Mat_SeqAIJ       *matA, *matB = NULL;
26*8be712e4SBarry Smith   Mat_MPIAIJ       *mpimat = NULL;
27*8be712e4SBarry Smith   MPI_Comm          comm;
28*8be712e4SBarry 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;
29*8be712e4SBarry Smith   PetscInt         *cpcol_gid, *cpcol_state, *lid_cprowID, *lid_gid, *cpcol_sel_gid, *icpcol_gid, *lid_state, *lid_parent_gid = NULL, nrm_tot = 0;
30*8be712e4SBarry Smith   PetscBool        *lid_removed;
31*8be712e4SBarry Smith   PetscBool         isMPI, isAIJ, isOK;
32*8be712e4SBarry Smith   const PetscInt   *perm_ix;
33*8be712e4SBarry Smith   const PetscInt    nloc = Gmat->rmap->n;
34*8be712e4SBarry Smith   PetscCoarsenData *agg_lists;
35*8be712e4SBarry Smith   PetscLayout       layout;
36*8be712e4SBarry Smith   PetscSF           sf;
37*8be712e4SBarry Smith   IS                info_is;
38*8be712e4SBarry Smith 
39*8be712e4SBarry Smith   PetscFunctionBegin;
40*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
41*8be712e4SBarry Smith   PetscCall(ISCreate(comm, &info_is));
42*8be712e4SBarry Smith   PetscCall(PetscInfo(info_is, "mis: nloc = %d\n", (int)nloc));
43*8be712e4SBarry Smith   /* get submatrices */
44*8be712e4SBarry Smith   PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPI));
45*8be712e4SBarry Smith   if (isMPI) {
46*8be712e4SBarry Smith     mpimat = (Mat_MPIAIJ *)Gmat->data;
47*8be712e4SBarry Smith     matA   = (Mat_SeqAIJ *)mpimat->A->data;
48*8be712e4SBarry Smith     matB   = (Mat_SeqAIJ *)mpimat->B->data;
49*8be712e4SBarry Smith     /* force compressed storage of B */
50*8be712e4SBarry Smith     PetscCall(MatCheckCompressedRow(mpimat->B, matB->nonzerorowcnt, &matB->compressedrow, matB->i, Gmat->rmap->n, -1.0));
51*8be712e4SBarry Smith   } else {
52*8be712e4SBarry Smith     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATSEQAIJ, &isAIJ));
53*8be712e4SBarry Smith     PetscCheck(isAIJ, comm, PETSC_ERR_PLIB, "Require AIJ matrix.");
54*8be712e4SBarry Smith     matA = (Mat_SeqAIJ *)Gmat->data;
55*8be712e4SBarry Smith   }
56*8be712e4SBarry Smith   PetscCall(MatGetOwnershipRange(Gmat, &my0, &Iend));
57*8be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_gid)); /* explicit array needed */
58*8be712e4SBarry Smith   if (mpimat) {
59*8be712e4SBarry Smith     for (kk = 0, gid = my0; kk < nloc; kk++, gid++) lid_gid[kk] = gid;
60*8be712e4SBarry Smith     PetscCall(VecGetLocalSize(mpimat->lvec, &num_fine_ghosts));
61*8be712e4SBarry Smith     PetscCall(PetscMalloc1(num_fine_ghosts, &cpcol_gid));
62*8be712e4SBarry Smith     PetscCall(PetscMalloc1(num_fine_ghosts, &cpcol_state));
63*8be712e4SBarry Smith     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)Gmat), &sf));
64*8be712e4SBarry Smith     PetscCall(MatGetLayouts(Gmat, &layout, NULL));
65*8be712e4SBarry Smith     PetscCall(PetscSFSetGraphLayout(sf, layout, num_fine_ghosts, NULL, PETSC_COPY_VALUES, mpimat->garray));
66*8be712e4SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
67*8be712e4SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
68*8be712e4SBarry Smith     for (kk = 0; kk < num_fine_ghosts; kk++) cpcol_state[kk] = MIS_NOT_DONE;
69*8be712e4SBarry Smith   } else num_fine_ghosts = 0;
70*8be712e4SBarry Smith 
71*8be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_cprowID));
72*8be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_removed)); /* explicit array needed */
73*8be712e4SBarry Smith   if (strict_aggs) PetscCall(PetscMalloc1(nloc, &lid_parent_gid));
74*8be712e4SBarry Smith   PetscCall(PetscMalloc1(nloc, &lid_state));
75*8be712e4SBarry Smith 
76*8be712e4SBarry Smith   /* has ghost nodes for !strict and uses local indexing (yuck) */
77*8be712e4SBarry Smith   PetscCall(PetscCDCreate(strict_aggs ? nloc : num_fine_ghosts + nloc, &agg_lists));
78*8be712e4SBarry Smith   if (a_locals_llist) *a_locals_llist = agg_lists;
79*8be712e4SBarry Smith 
80*8be712e4SBarry Smith   /* need an inverse map - locals */
81*8be712e4SBarry Smith   for (kk = 0; kk < nloc; kk++) {
82*8be712e4SBarry Smith     lid_cprowID[kk] = -1;
83*8be712e4SBarry Smith     lid_removed[kk] = PETSC_FALSE;
84*8be712e4SBarry Smith     if (strict_aggs) lid_parent_gid[kk] = -1.0;
85*8be712e4SBarry Smith     lid_state[kk] = MIS_NOT_DONE;
86*8be712e4SBarry Smith   }
87*8be712e4SBarry Smith   /* set index into cmpressed row 'lid_cprowID' */
88*8be712e4SBarry Smith   if (matB) {
89*8be712e4SBarry Smith     for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
90*8be712e4SBarry Smith       lid = matB->compressedrow.rindex[ix];
91*8be712e4SBarry Smith       if (lid >= 0) lid_cprowID[lid] = ix;
92*8be712e4SBarry Smith     }
93*8be712e4SBarry Smith   }
94*8be712e4SBarry Smith   /* MIS */
95*8be712e4SBarry Smith   nremoved = nDone = 0;
96*8be712e4SBarry Smith 
97*8be712e4SBarry Smith   PetscCall(ISGetIndices(perm, &perm_ix));
98*8be712e4SBarry Smith   while (nDone < nloc || PETSC_TRUE) { /* asynchronous not implemented */
99*8be712e4SBarry Smith     /* check all vertices */
100*8be712e4SBarry Smith     for (kk = 0; kk < nloc; kk++) {
101*8be712e4SBarry Smith       lid   = perm_ix[kk];
102*8be712e4SBarry Smith       state = lid_state[lid];
103*8be712e4SBarry Smith       if (lid_removed[lid]) continue;
104*8be712e4SBarry Smith       if (state == MIS_NOT_DONE) {
105*8be712e4SBarry Smith         /* parallel test, delete if selected ghost */
106*8be712e4SBarry Smith         isOK = PETSC_TRUE;
107*8be712e4SBarry Smith         if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
108*8be712e4SBarry Smith           ii  = matB->compressedrow.i;
109*8be712e4SBarry Smith           n   = ii[ix + 1] - ii[ix];
110*8be712e4SBarry Smith           idx = matB->j + ii[ix];
111*8be712e4SBarry Smith           for (j = 0; j < n; j++) {
112*8be712e4SBarry Smith             cpid   = idx[j]; /* compressed row ID in B mat */
113*8be712e4SBarry Smith             gid    = cpcol_gid[cpid];
114*8be712e4SBarry Smith             statej = cpcol_state[cpid];
115*8be712e4SBarry Smith             PetscCheck(!MIS_IS_SELECTED(statej), PETSC_COMM_SELF, PETSC_ERR_SUP, "selected ghost: %d", (int)gid);
116*8be712e4SBarry Smith             if (statej == MIS_NOT_DONE && gid >= Iend) { /* should be (pe>rank), use gid as pe proxy */
117*8be712e4SBarry Smith               isOK = PETSC_FALSE;                        /* can not delete */
118*8be712e4SBarry Smith               break;
119*8be712e4SBarry Smith             }
120*8be712e4SBarry Smith           }
121*8be712e4SBarry Smith         }           /* parallel test */
122*8be712e4SBarry Smith         if (isOK) { /* select or remove this vertex */
123*8be712e4SBarry Smith           nDone++;
124*8be712e4SBarry Smith           /* check for singleton */
125*8be712e4SBarry Smith           ii = matA->i;
126*8be712e4SBarry Smith           n  = ii[lid + 1] - ii[lid];
127*8be712e4SBarry Smith           if (n < 2) {
128*8be712e4SBarry Smith             /* if I have any ghost adj then not a sing */
129*8be712e4SBarry Smith             ix = lid_cprowID[lid];
130*8be712e4SBarry Smith             if (ix == -1 || !(matB->compressedrow.i[ix + 1] - matB->compressedrow.i[ix])) {
131*8be712e4SBarry Smith               nremoved++;
132*8be712e4SBarry Smith               nrm_tot++;
133*8be712e4SBarry Smith               lid_removed[lid] = PETSC_TRUE;
134*8be712e4SBarry Smith               continue;
135*8be712e4SBarry Smith               // lid_state[lidj] = MIS_REMOVED; /* add singleton to MIS (can cause low rank with elasticity on fine grid) */
136*8be712e4SBarry Smith             }
137*8be712e4SBarry Smith           }
138*8be712e4SBarry Smith           /* SELECTED state encoded with global index */
139*8be712e4SBarry Smith           lid_state[lid] = lid + my0;
140*8be712e4SBarry Smith           nselected++;
141*8be712e4SBarry Smith           if (strict_aggs) {
142*8be712e4SBarry Smith             PetscCall(PetscCDAppendID(agg_lists, lid, lid + my0));
143*8be712e4SBarry Smith           } else {
144*8be712e4SBarry Smith             PetscCall(PetscCDAppendID(agg_lists, lid, lid));
145*8be712e4SBarry Smith           }
146*8be712e4SBarry Smith           /* delete local adj */
147*8be712e4SBarry Smith           idx = matA->j + ii[lid];
148*8be712e4SBarry Smith           for (j = 0; j < n; j++) {
149*8be712e4SBarry Smith             lidj   = idx[j];
150*8be712e4SBarry Smith             statej = lid_state[lidj];
151*8be712e4SBarry Smith             if (statej == MIS_NOT_DONE) {
152*8be712e4SBarry Smith               nDone++;
153*8be712e4SBarry Smith               if (strict_aggs) {
154*8be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lid, lidj + my0));
155*8be712e4SBarry Smith               } else {
156*8be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lid, lidj));
157*8be712e4SBarry Smith               }
158*8be712e4SBarry Smith               lid_state[lidj] = MIS_DELETED; /* delete this */
159*8be712e4SBarry Smith             }
160*8be712e4SBarry Smith           }
161*8be712e4SBarry Smith           /* delete ghost adj of lid - deleted ghost done later for strict_aggs */
162*8be712e4SBarry Smith           if (!strict_aggs) {
163*8be712e4SBarry Smith             if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
164*8be712e4SBarry Smith               ii  = matB->compressedrow.i;
165*8be712e4SBarry Smith               n   = ii[ix + 1] - ii[ix];
166*8be712e4SBarry Smith               idx = matB->j + ii[ix];
167*8be712e4SBarry Smith               for (j = 0; j < n; j++) {
168*8be712e4SBarry Smith                 cpid   = idx[j]; /* compressed row ID in B mat */
169*8be712e4SBarry Smith                 statej = cpcol_state[cpid];
170*8be712e4SBarry Smith                 if (statej == MIS_NOT_DONE) PetscCall(PetscCDAppendID(agg_lists, lid, nloc + cpid));
171*8be712e4SBarry Smith               }
172*8be712e4SBarry Smith             }
173*8be712e4SBarry Smith           }
174*8be712e4SBarry Smith         } /* selected */
175*8be712e4SBarry Smith       }   /* not done vertex */
176*8be712e4SBarry Smith     }     /* vertex loop */
177*8be712e4SBarry Smith 
178*8be712e4SBarry Smith     /* update ghost states and count todos */
179*8be712e4SBarry Smith     if (mpimat) {
180*8be712e4SBarry Smith       /* scatter states, check for done */
181*8be712e4SBarry Smith       PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
182*8be712e4SBarry Smith       PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
183*8be712e4SBarry Smith       ii = matB->compressedrow.i;
184*8be712e4SBarry Smith       for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
185*8be712e4SBarry Smith         lid   = matB->compressedrow.rindex[ix]; /* local boundary node */
186*8be712e4SBarry Smith         state = lid_state[lid];
187*8be712e4SBarry Smith         if (state == MIS_NOT_DONE) {
188*8be712e4SBarry Smith           /* look at ghosts */
189*8be712e4SBarry Smith           n   = ii[ix + 1] - ii[ix];
190*8be712e4SBarry Smith           idx = matB->j + ii[ix];
191*8be712e4SBarry Smith           for (j = 0; j < n; j++) {
192*8be712e4SBarry Smith             cpid   = idx[j]; /* compressed row ID in B mat */
193*8be712e4SBarry Smith             statej = cpcol_state[cpid];
194*8be712e4SBarry Smith             if (MIS_IS_SELECTED(statej)) { /* lid is now deleted, do it */
195*8be712e4SBarry Smith               nDone++;
196*8be712e4SBarry Smith               lid_state[lid] = MIS_DELETED; /* delete this */
197*8be712e4SBarry Smith               if (!strict_aggs) {
198*8be712e4SBarry Smith                 lidj = nloc + cpid;
199*8be712e4SBarry Smith                 PetscCall(PetscCDAppendID(agg_lists, lidj, lid));
200*8be712e4SBarry Smith               } else {
201*8be712e4SBarry Smith                 sgid                = cpcol_gid[cpid];
202*8be712e4SBarry Smith                 lid_parent_gid[lid] = sgid; /* keep track of proc that I belong to */
203*8be712e4SBarry Smith               }
204*8be712e4SBarry Smith               break;
205*8be712e4SBarry Smith             }
206*8be712e4SBarry Smith           }
207*8be712e4SBarry Smith         }
208*8be712e4SBarry Smith       }
209*8be712e4SBarry Smith       /* all done? */
210*8be712e4SBarry Smith       t1 = nloc - nDone;
211*8be712e4SBarry Smith       PetscCall(MPIU_Allreduce(&t1, &t2, 1, MPIU_INT, MPI_SUM, comm)); /* synchronous version */
212*8be712e4SBarry Smith       if (!t2) break;
213*8be712e4SBarry Smith     } else break; /* all done */
214*8be712e4SBarry Smith   }               /* outer parallel MIS loop */
215*8be712e4SBarry Smith   PetscCall(ISRestoreIndices(perm, &perm_ix));
216*8be712e4SBarry Smith   PetscCall(PetscInfo(info_is, "\t removed %" PetscInt_FMT " of %" PetscInt_FMT " vertices.  %" PetscInt_FMT " selected.\n", nremoved, nloc, nselected));
217*8be712e4SBarry Smith 
218*8be712e4SBarry Smith   /* tell adj who my lid_parent_gid vertices belong to - fill in agg_lists selected ghost lists */
219*8be712e4SBarry Smith   if (strict_aggs && matB) {
220*8be712e4SBarry Smith     /* need to copy this to free buffer -- should do this globally */
221*8be712e4SBarry Smith     PetscCall(PetscMalloc1(num_fine_ghosts, &cpcol_sel_gid));
222*8be712e4SBarry Smith     PetscCall(PetscMalloc1(num_fine_ghosts, &icpcol_gid));
223*8be712e4SBarry Smith     for (cpid = 0; cpid < num_fine_ghosts; cpid++) icpcol_gid[cpid] = cpcol_gid[cpid];
224*8be712e4SBarry Smith 
225*8be712e4SBarry Smith     /* get proc of deleted ghost */
226*8be712e4SBarry Smith     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
227*8be712e4SBarry Smith     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
228*8be712e4SBarry Smith     for (cpid = 0; cpid < num_fine_ghosts; cpid++) {
229*8be712e4SBarry Smith       sgid = cpcol_sel_gid[cpid];
230*8be712e4SBarry Smith       gid  = icpcol_gid[cpid];
231*8be712e4SBarry Smith       if (sgid >= my0 && sgid < Iend) { /* I own this deleted */
232*8be712e4SBarry Smith         slid = sgid - my0;
233*8be712e4SBarry Smith         PetscCall(PetscCDAppendID(agg_lists, slid, gid));
234*8be712e4SBarry Smith       }
235*8be712e4SBarry Smith     }
236*8be712e4SBarry Smith     PetscCall(PetscFree(icpcol_gid));
237*8be712e4SBarry Smith     PetscCall(PetscFree(cpcol_sel_gid));
238*8be712e4SBarry Smith   }
239*8be712e4SBarry Smith   if (mpimat) {
240*8be712e4SBarry Smith     PetscCall(PetscSFDestroy(&sf));
241*8be712e4SBarry Smith     PetscCall(PetscFree(cpcol_gid));
242*8be712e4SBarry Smith     PetscCall(PetscFree(cpcol_state));
243*8be712e4SBarry Smith   }
244*8be712e4SBarry Smith   PetscCall(PetscFree(lid_cprowID));
245*8be712e4SBarry Smith   PetscCall(PetscFree(lid_gid));
246*8be712e4SBarry Smith   PetscCall(PetscFree(lid_removed));
247*8be712e4SBarry Smith   if (strict_aggs) PetscCall(PetscFree(lid_parent_gid));
248*8be712e4SBarry Smith   PetscCall(PetscFree(lid_state));
249*8be712e4SBarry Smith   if (strict_aggs) {
250*8be712e4SBarry Smith     // check sizes -- all vertices must get in graph
251*8be712e4SBarry Smith     PetscInt aa[2] = {0, nrm_tot}, bb[2], MM;
252*8be712e4SBarry Smith     PetscCall(MatGetSize(Gmat, &MM, NULL));
253*8be712e4SBarry Smith     // check sizes -- all vertices must get in graph
254*8be712e4SBarry Smith     PetscCall(PetscCDCount(agg_lists, &aa[0]));
255*8be712e4SBarry Smith     PetscCall(MPIU_Allreduce(aa, bb, 2, MPIU_INT, MPI_SUM, comm));
256*8be712e4SBarry 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]));
257*8be712e4SBarry Smith     PetscCheck(MM >= bb[0], comm, PETSC_ERR_PLIB, "Sum of aggs too big");
258*8be712e4SBarry Smith   }
259*8be712e4SBarry Smith   PetscCall(ISDestroy(&info_is));
260*8be712e4SBarry Smith 
261*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
262*8be712e4SBarry Smith }
263*8be712e4SBarry Smith 
264*8be712e4SBarry Smith /*
265*8be712e4SBarry Smith    MIS coarsen, simple greedy.
266*8be712e4SBarry Smith */
267*8be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MIS(MatCoarsen coarse)
268*8be712e4SBarry Smith {
269*8be712e4SBarry Smith   Mat mat = coarse->graph;
270*8be712e4SBarry Smith 
271*8be712e4SBarry Smith   PetscFunctionBegin;
272*8be712e4SBarry Smith   if (!coarse->perm) {
273*8be712e4SBarry Smith     IS       perm;
274*8be712e4SBarry Smith     PetscInt n, m;
275*8be712e4SBarry Smith     MPI_Comm comm;
276*8be712e4SBarry Smith 
277*8be712e4SBarry Smith     PetscCall(PetscObjectGetComm((PetscObject)mat, &comm));
278*8be712e4SBarry Smith     PetscCall(MatGetLocalSize(mat, &m, &n));
279*8be712e4SBarry Smith     PetscCall(ISCreateStride(comm, m, 0, 1, &perm));
280*8be712e4SBarry Smith     PetscCall(MatCoarsenApply_MIS_private(perm, mat, coarse->strict_aggs, &coarse->agg_lists));
281*8be712e4SBarry Smith     PetscCall(ISDestroy(&perm));
282*8be712e4SBarry Smith   } else {
283*8be712e4SBarry Smith     PetscCall(MatCoarsenApply_MIS_private(coarse->perm, mat, coarse->strict_aggs, &coarse->agg_lists));
284*8be712e4SBarry Smith   }
285*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
286*8be712e4SBarry Smith }
287*8be712e4SBarry Smith 
288*8be712e4SBarry Smith static PetscErrorCode MatCoarsenView_MIS(MatCoarsen coarse, PetscViewer viewer)
289*8be712e4SBarry Smith {
290*8be712e4SBarry Smith   PetscMPIInt rank;
291*8be712e4SBarry Smith   PetscBool   iascii;
292*8be712e4SBarry Smith 
293*8be712e4SBarry Smith   PetscFunctionBegin;
294*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank));
295*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
296*8be712e4SBarry Smith   if (iascii) {
297*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
298*8be712e4SBarry Smith     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] MIS aggregator\n", rank));
299*8be712e4SBarry Smith     if (!rank) {
300*8be712e4SBarry Smith       PetscCDIntNd *pos, *pos2;
301*8be712e4SBarry Smith       for (PetscInt kk = 0; kk < coarse->agg_lists->size; kk++) {
302*8be712e4SBarry Smith         PetscCall(PetscCDGetHeadPos(coarse->agg_lists, kk, &pos));
303*8be712e4SBarry Smith         if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected %d: ", (int)kk));
304*8be712e4SBarry Smith         while (pos) {
305*8be712e4SBarry Smith           PetscInt gid1;
306*8be712e4SBarry Smith           PetscCall(PetscCDIntNdGetID(pos, &gid1));
307*8be712e4SBarry Smith           PetscCall(PetscCDGetNextPos(coarse->agg_lists, kk, &pos));
308*8be712e4SBarry Smith           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %d ", (int)gid1));
309*8be712e4SBarry Smith         }
310*8be712e4SBarry Smith         if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
311*8be712e4SBarry Smith       }
312*8be712e4SBarry Smith     }
313*8be712e4SBarry Smith     PetscCall(PetscViewerFlush(viewer));
314*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
315*8be712e4SBarry Smith   }
316*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
317*8be712e4SBarry Smith }
318*8be712e4SBarry Smith 
319*8be712e4SBarry Smith /*MC
320*8be712e4SBarry Smith    MATCOARSENMIS - Creates a coarsening with a maximal independent set (MIS) algorithm
321*8be712e4SBarry Smith 
322*8be712e4SBarry Smith    Collective
323*8be712e4SBarry Smith 
324*8be712e4SBarry Smith    Input Parameter:
325*8be712e4SBarry Smith .  coarse - the coarsen context
326*8be712e4SBarry Smith 
327*8be712e4SBarry Smith    Level: beginner
328*8be712e4SBarry Smith 
329*8be712e4SBarry Smith .seealso: `MatCoarsen`, `MatCoarsenApply()`, `MatCoarsenGetData()`, `MatCoarsenSetType()`, `MatCoarsenType`
330*8be712e4SBarry Smith M*/
331*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_MIS(MatCoarsen coarse)
332*8be712e4SBarry Smith {
333*8be712e4SBarry Smith   PetscFunctionBegin;
334*8be712e4SBarry Smith   coarse->ops->apply = MatCoarsenApply_MIS;
335*8be712e4SBarry Smith   coarse->ops->view  = MatCoarsenView_MIS;
336*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
337*8be712e4SBarry Smith }
338