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 >= 0)
108be712e4SBarry Smith
118be712e4SBarry Smith /* edge for priority queue */
128be712e4SBarry Smith typedef struct edge_tag {
138be712e4SBarry Smith PetscReal weight;
148be712e4SBarry Smith PetscInt lid0, gid1, cpid1;
158be712e4SBarry Smith } Edge;
168be712e4SBarry Smith
PetscCoarsenDataView_private(PetscCoarsenData * agg_lists,PetscViewer viewer)178be712e4SBarry Smith static PetscErrorCode PetscCoarsenDataView_private(PetscCoarsenData *agg_lists, PetscViewer viewer)
188be712e4SBarry Smith {
198be712e4SBarry Smith PetscCDIntNd *pos, *pos2;
208be712e4SBarry Smith
218be712e4SBarry Smith PetscFunctionBegin;
228be712e4SBarry Smith for (PetscInt kk = 0; kk < agg_lists->size; kk++) {
238be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_lists, kk, &pos));
24835f2295SStefano Zampini if ((pos2 = pos)) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "selected local %" PetscInt_FMT ": ", kk));
258be712e4SBarry Smith while (pos) {
268be712e4SBarry Smith PetscInt gid1;
278be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gid1));
288be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_lists, kk, &pos));
29835f2295SStefano Zampini PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT " ", gid1));
308be712e4SBarry Smith }
318be712e4SBarry Smith if (pos2) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n"));
328be712e4SBarry Smith }
338be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
348be712e4SBarry Smith }
358be712e4SBarry Smith
368be712e4SBarry Smith /*
378be712e4SBarry Smith MatCoarsenApply_MISK_private - parallel heavy edge matching
388be712e4SBarry Smith
398be712e4SBarry Smith Input Parameter:
408be712e4SBarry Smith . perm - permutation
418be712e4SBarry Smith . Gmat - global matrix of graph (data not defined)
428be712e4SBarry Smith
438be712e4SBarry Smith Output Parameter:
448be712e4SBarry Smith . a_locals_llist - array of list of local nodes rooted at local node
458be712e4SBarry Smith */
MatCoarsenApply_MISK_private(IS perm,const PetscInt misk,Mat Gmat,PetscCoarsenData ** a_locals_llist)468be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MISK_private(IS perm, const PetscInt misk, Mat Gmat, PetscCoarsenData **a_locals_llist)
478be712e4SBarry Smith {
488be712e4SBarry Smith PetscBool isMPI;
498be712e4SBarry Smith MPI_Comm comm;
508be712e4SBarry Smith PetscMPIInt rank, size;
518be712e4SBarry Smith Mat cMat, Prols[5], Rtot;
528be712e4SBarry Smith PetscScalar one = 1;
538be712e4SBarry Smith
548be712e4SBarry Smith PetscFunctionBegin;
558be712e4SBarry Smith PetscValidHeaderSpecific(perm, IS_CLASSID, 1);
568be712e4SBarry Smith PetscValidHeaderSpecific(Gmat, MAT_CLASSID, 3);
578be712e4SBarry Smith PetscAssertPointer(a_locals_llist, 4);
58835f2295SStefano Zampini PetscCheck(misk < 5 && misk > 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "too many/few levels: %" PetscInt_FMT, misk);
598be712e4SBarry Smith PetscCall(PetscObjectBaseTypeCompare((PetscObject)Gmat, MATMPIAIJ, &isMPI));
608be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)Gmat, &comm));
618be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
628be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size));
63835f2295SStefano Zampini PetscCall(PetscInfo(Gmat, "misk %" PetscInt_FMT "\n", misk));
648be712e4SBarry Smith /* make a copy of the graph, this gets destroyed in iterates */
658be712e4SBarry Smith if (misk > 1) PetscCall(MatDuplicate(Gmat, MAT_COPY_VALUES, &cMat));
668be712e4SBarry Smith else cMat = Gmat;
678be712e4SBarry Smith for (PetscInt iterIdx = 0; iterIdx < misk; iterIdx++) {
688be712e4SBarry Smith Mat_SeqAIJ *matA, *matB = NULL;
698be712e4SBarry Smith Mat_MPIAIJ *mpimat = NULL;
708be712e4SBarry Smith const PetscInt *perm_ix;
718be712e4SBarry Smith const PetscInt nloc_inner = cMat->rmap->n;
728be712e4SBarry Smith PetscCoarsenData *agg_lists;
738be712e4SBarry Smith PetscInt *cpcol_gid = NULL, *cpcol_state, *lid_cprowID, *lid_state, *lid_parent_gid = NULL;
748be712e4SBarry Smith PetscInt num_fine_ghosts, kk, n, ix, j, *idx, *ai, Iend, my0, nremoved, gid, cpid, lidj, sgid, t1, t2, slid, nDone, nselected = 0, state;
758be712e4SBarry Smith PetscBool *lid_removed, isOK;
768be712e4SBarry Smith PetscLayout layout;
778be712e4SBarry Smith PetscSF sf;
788be712e4SBarry Smith
798be712e4SBarry Smith if (isMPI) {
808be712e4SBarry Smith mpimat = (Mat_MPIAIJ *)cMat->data;
818be712e4SBarry Smith matA = (Mat_SeqAIJ *)mpimat->A->data;
828be712e4SBarry Smith matB = (Mat_SeqAIJ *)mpimat->B->data;
838be712e4SBarry Smith /* force compressed storage of B */
848be712e4SBarry Smith PetscCall(MatCheckCompressedRow(mpimat->B, matB->nonzerorowcnt, &matB->compressedrow, matB->i, cMat->rmap->n, -1.0));
858be712e4SBarry Smith } else {
868be712e4SBarry Smith PetscBool isAIJ;
87c376f201SBarry Smith
88c376f201SBarry Smith matA = (Mat_SeqAIJ *)cMat->data;
898be712e4SBarry Smith PetscCall(PetscObjectBaseTypeCompare((PetscObject)cMat, MATSEQAIJ, &isAIJ));
908be712e4SBarry Smith PetscCheck(isAIJ, PETSC_COMM_SELF, PETSC_ERR_USER, "Require AIJ matrix.");
918be712e4SBarry Smith }
928be712e4SBarry Smith PetscCall(MatGetOwnershipRange(cMat, &my0, &Iend));
93c376f201SBarry Smith if (isMPI) {
948be712e4SBarry Smith PetscInt *lid_gid;
95c376f201SBarry Smith
968be712e4SBarry Smith PetscCall(PetscMalloc1(nloc_inner, &lid_gid)); /* explicit array needed */
978be712e4SBarry Smith for (kk = 0, gid = my0; kk < nloc_inner; kk++, gid++) lid_gid[kk] = gid;
988be712e4SBarry Smith PetscCall(VecGetLocalSize(mpimat->lvec, &num_fine_ghosts));
99c376f201SBarry Smith PetscCall(PetscMalloc2(num_fine_ghosts, &cpcol_gid, num_fine_ghosts, &cpcol_state));
1008be712e4SBarry Smith PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)cMat), &sf));
1018be712e4SBarry Smith PetscCall(MatGetLayouts(cMat, &layout, NULL));
1028be712e4SBarry Smith PetscCall(PetscSFSetGraphLayout(sf, layout, num_fine_ghosts, NULL, PETSC_COPY_VALUES, mpimat->garray));
1038be712e4SBarry Smith PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
1048be712e4SBarry Smith PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_gid, cpcol_gid, MPI_REPLACE));
1058be712e4SBarry Smith for (kk = 0; kk < num_fine_ghosts; kk++) cpcol_state[kk] = MIS_NOT_DONE;
1068be712e4SBarry Smith PetscCall(PetscFree(lid_gid));
1078be712e4SBarry Smith } else num_fine_ghosts = 0;
1088be712e4SBarry Smith
109c376f201SBarry Smith PetscCall(PetscMalloc4(nloc_inner, &lid_cprowID, nloc_inner, &lid_removed, nloc_inner, &lid_parent_gid, nloc_inner, &lid_state));
1108be712e4SBarry Smith PetscCall(PetscCDCreate(nloc_inner, &agg_lists));
1118be712e4SBarry Smith /* need an inverse map - locals */
1128be712e4SBarry Smith for (kk = 0; kk < nloc_inner; kk++) {
1138be712e4SBarry Smith lid_cprowID[kk] = -1;
1148be712e4SBarry Smith lid_removed[kk] = PETSC_FALSE;
1158be712e4SBarry Smith lid_parent_gid[kk] = -1.0;
1168be712e4SBarry Smith lid_state[kk] = MIS_NOT_DONE;
1178be712e4SBarry Smith }
1188be712e4SBarry Smith /* set index into cmpressed row 'lid_cprowID' */
1198be712e4SBarry Smith if (matB) {
1208be712e4SBarry Smith for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
1218be712e4SBarry Smith const PetscInt lid = matB->compressedrow.rindex[ix];
1228be712e4SBarry Smith if (lid >= 0) lid_cprowID[lid] = ix;
1238be712e4SBarry Smith }
1248be712e4SBarry Smith }
1258be712e4SBarry Smith /* MIS */
1268be712e4SBarry Smith nremoved = nDone = 0;
1278be712e4SBarry Smith if (!iterIdx) PetscCall(ISGetIndices(perm, &perm_ix)); // use permutation on first MIS
1288be712e4SBarry Smith else perm_ix = NULL;
1298be712e4SBarry Smith while (nDone < nloc_inner || PETSC_TRUE) { /* asynchronous not implemented */
1308be712e4SBarry Smith /* check all vertices */
1318be712e4SBarry Smith for (kk = 0; kk < nloc_inner; kk++) {
1328be712e4SBarry Smith const PetscInt lid = perm_ix ? perm_ix[kk] : kk;
1338be712e4SBarry Smith state = lid_state[lid];
1348be712e4SBarry Smith if (iterIdx == 0 && lid_removed[lid]) continue;
1358be712e4SBarry Smith if (state == MIS_NOT_DONE) {
1368be712e4SBarry Smith /* parallel test, delete if selected ghost */
1378be712e4SBarry Smith isOK = PETSC_TRUE;
1388be712e4SBarry Smith /* parallel test */
1398be712e4SBarry Smith if ((ix = lid_cprowID[lid]) != -1) { /* if I have any ghost neighbors */
1408be712e4SBarry Smith ai = matB->compressedrow.i;
1418be712e4SBarry Smith n = ai[ix + 1] - ai[ix];
1428be712e4SBarry Smith idx = matB->j + ai[ix];
1438be712e4SBarry Smith for (j = 0; j < n; j++) {
1448be712e4SBarry Smith cpid = idx[j]; /* compressed row ID in B mat */
1458be712e4SBarry Smith gid = cpcol_gid[cpid];
1468be712e4SBarry Smith if (cpcol_state[cpid] == MIS_NOT_DONE && gid >= Iend) { /* or pe>rank */
1478be712e4SBarry Smith isOK = PETSC_FALSE; /* can not delete */
1488be712e4SBarry Smith break;
1498be712e4SBarry Smith }
1508be712e4SBarry Smith }
1518be712e4SBarry Smith }
1528be712e4SBarry Smith if (isOK) { /* select or remove this vertex if it is a true singleton like a BC */
1538be712e4SBarry Smith nDone++;
1548be712e4SBarry Smith /* check for singleton */
1558be712e4SBarry Smith ai = matA->i;
1568be712e4SBarry Smith n = ai[lid + 1] - ai[lid];
1578be712e4SBarry Smith if (n < 2) {
1588be712e4SBarry Smith /* if I have any ghost adj then not a singleton */
1598be712e4SBarry Smith ix = lid_cprowID[lid];
1608be712e4SBarry Smith if (ix == -1 || !(matB->compressedrow.i[ix + 1] - matB->compressedrow.i[ix])) {
1618be712e4SBarry Smith if (iterIdx == 0) {
1628be712e4SBarry Smith lid_removed[lid] = PETSC_TRUE;
1638be712e4SBarry Smith nremoved++; // let it get selected
1648be712e4SBarry Smith }
1658be712e4SBarry Smith // PetscCall(PetscCDAppendID(agg_lists, lid, lid + my0));
1668be712e4SBarry Smith // lid_state[lid] = nselected; // >= 0 is selected, cache for ordering coarse grid
1678be712e4SBarry Smith /* should select this because it is technically in the MIS but lets not */
1688be712e4SBarry Smith continue; /* one local adj (me) and no ghost - singleton */
1698be712e4SBarry Smith }
1708be712e4SBarry Smith }
1718be712e4SBarry Smith /* SELECTED state encoded with global index */
1728be712e4SBarry Smith lid_state[lid] = nselected; // >= 0 is selected, cache for ordering coarse grid
1738be712e4SBarry Smith nselected++;
1748be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_lists, lid, lid + my0));
1758be712e4SBarry Smith /* delete local adj */
1768be712e4SBarry Smith idx = matA->j + ai[lid];
1778be712e4SBarry Smith for (j = 0; j < n; j++) {
1788be712e4SBarry Smith lidj = idx[j];
1798be712e4SBarry Smith if (lid_state[lidj] == MIS_NOT_DONE) {
1808be712e4SBarry Smith nDone++;
1818be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_lists, lid, lidj + my0));
1828be712e4SBarry Smith lid_state[lidj] = MIS_DELETED; /* delete this */
1838be712e4SBarry Smith }
1848be712e4SBarry Smith }
1858be712e4SBarry Smith } /* selected */
1868be712e4SBarry Smith } /* not done vertex */
1878be712e4SBarry Smith } /* vertex loop */
1888be712e4SBarry Smith
1898be712e4SBarry Smith /* update ghost states and count todos */
190c376f201SBarry Smith if (isMPI) {
1918be712e4SBarry Smith /* scatter states, check for done */
1928be712e4SBarry Smith PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
1938be712e4SBarry Smith PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_state, cpcol_state, MPI_REPLACE));
1948be712e4SBarry Smith ai = matB->compressedrow.i;
1958be712e4SBarry Smith for (ix = 0; ix < matB->compressedrow.nrows; ix++) {
1966497c311SBarry Smith const PetscInt lidj = matB->compressedrow.rindex[ix]; /* local boundary node */
1978be712e4SBarry Smith state = lid_state[lidj];
1988be712e4SBarry Smith if (state == MIS_NOT_DONE) {
1998be712e4SBarry Smith /* look at ghosts */
2008be712e4SBarry Smith n = ai[ix + 1] - ai[ix];
2018be712e4SBarry Smith idx = matB->j + ai[ix];
2028be712e4SBarry Smith for (j = 0; j < n; j++) {
2038be712e4SBarry Smith cpid = idx[j]; /* compressed row ID in B mat */
2048be712e4SBarry Smith if (MIS_IS_SELECTED(cpcol_state[cpid])) { /* lid is now deleted by ghost */
2058be712e4SBarry Smith nDone++;
2068be712e4SBarry Smith lid_state[lidj] = MIS_DELETED; /* delete this */
2078be712e4SBarry Smith sgid = cpcol_gid[cpid];
2088be712e4SBarry Smith lid_parent_gid[lidj] = sgid; /* keep track of proc that I belong to */
2098be712e4SBarry Smith break;
2108be712e4SBarry Smith }
2118be712e4SBarry Smith }
2128be712e4SBarry Smith }
2138be712e4SBarry Smith }
2148be712e4SBarry Smith /* all done? */
2158be712e4SBarry Smith t1 = nloc_inner - nDone;
216462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&t1, &t2, 1, MPIU_INT, MPI_SUM, comm)); /* synchronous version */
2178be712e4SBarry Smith if (!t2) break;
2188be712e4SBarry Smith } else break; /* no mpi - all done */
2198be712e4SBarry Smith } /* outer parallel MIS loop */
2208be712e4SBarry Smith if (!iterIdx) PetscCall(ISRestoreIndices(perm, &perm_ix));
2218be712e4SBarry Smith PetscCall(PetscInfo(Gmat, "\t removed %" PetscInt_FMT " of %" PetscInt_FMT " vertices. %" PetscInt_FMT " selected.\n", nremoved, nloc_inner, nselected));
2228be712e4SBarry Smith
2238be712e4SBarry Smith /* tell adj who my lid_parent_gid vertices belong to - fill in agg_lists selected ghost lists */
2248be712e4SBarry Smith if (matB) {
2258be712e4SBarry Smith PetscInt *cpcol_sel_gid, *icpcol_gid;
226c376f201SBarry Smith
2278be712e4SBarry Smith /* need to copy this to free buffer -- should do this globally */
228c376f201SBarry Smith PetscCall(PetscMalloc2(num_fine_ghosts, &icpcol_gid, num_fine_ghosts, &cpcol_sel_gid));
2298be712e4SBarry Smith for (cpid = 0; cpid < num_fine_ghosts; cpid++) icpcol_gid[cpid] = cpcol_gid[cpid];
2308be712e4SBarry Smith /* get proc of deleted ghost */
2318be712e4SBarry Smith PetscCall(PetscSFBcastBegin(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
2328be712e4SBarry Smith PetscCall(PetscSFBcastEnd(sf, MPIU_INT, lid_parent_gid, cpcol_sel_gid, MPI_REPLACE));
2338be712e4SBarry Smith for (cpid = 0; cpid < num_fine_ghosts; cpid++) {
2348be712e4SBarry Smith sgid = cpcol_sel_gid[cpid];
2358be712e4SBarry Smith gid = icpcol_gid[cpid];
2368be712e4SBarry Smith if (sgid >= my0 && sgid < Iend) { /* I own this deleted */
2378be712e4SBarry Smith slid = sgid - my0;
2388be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_lists, slid, gid));
2398be712e4SBarry Smith }
2408be712e4SBarry Smith }
2418be712e4SBarry Smith // done - cleanup
242c376f201SBarry Smith PetscCall(PetscFree2(icpcol_gid, cpcol_sel_gid));
2438be712e4SBarry Smith PetscCall(PetscSFDestroy(&sf));
244c376f201SBarry Smith PetscCall(PetscFree2(cpcol_gid, cpcol_state));
2458be712e4SBarry Smith }
246c376f201SBarry Smith PetscCall(PetscFree4(lid_cprowID, lid_removed, lid_parent_gid, lid_state));
2478be712e4SBarry Smith
2488be712e4SBarry Smith /* MIS done - make projection matrix - P */
2498be712e4SBarry Smith MatType jtype;
2508be712e4SBarry Smith PetscCall(MatGetType(Gmat, &jtype));
2518be712e4SBarry Smith PetscCall(MatCreate(comm, &Prols[iterIdx]));
2528be712e4SBarry Smith PetscCall(MatSetType(Prols[iterIdx], jtype));
2538be712e4SBarry Smith PetscCall(MatSetSizes(Prols[iterIdx], nloc_inner, nselected, PETSC_DETERMINE, PETSC_DETERMINE));
2548be712e4SBarry Smith PetscCall(MatSeqAIJSetPreallocation(Prols[iterIdx], 1, NULL));
2558be712e4SBarry Smith PetscCall(MatMPIAIJSetPreallocation(Prols[iterIdx], 1, NULL, 1, NULL));
2568be712e4SBarry Smith {
2578be712e4SBarry Smith PetscCDIntNd *pos, *pos2;
2588be712e4SBarry Smith PetscInt colIndex, Iend, fgid;
259c376f201SBarry Smith
2608be712e4SBarry Smith PetscCall(MatGetOwnershipRangeColumn(Prols[iterIdx], &colIndex, &Iend));
2618be712e4SBarry Smith // TODO - order with permutation in lid_selected (reversed)
2628be712e4SBarry Smith for (PetscInt lid = 0; lid < agg_lists->size; lid++) {
2638be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_lists, lid, &pos));
2648be712e4SBarry Smith pos2 = pos;
2658be712e4SBarry Smith while (pos) {
2668be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &fgid));
2678be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_lists, lid, &pos));
2688be712e4SBarry Smith PetscCall(MatSetValues(Prols[iterIdx], 1, &fgid, 1, &colIndex, &one, INSERT_VALUES));
2698be712e4SBarry Smith }
2708be712e4SBarry Smith if (pos2) colIndex++;
2718be712e4SBarry Smith }
272835f2295SStefano Zampini PetscCheck(Iend == colIndex, PETSC_COMM_SELF, PETSC_ERR_SUP, "Iend!=colIndex: %" PetscInt_FMT " %" PetscInt_FMT, Iend, colIndex);
2738be712e4SBarry Smith }
2748be712e4SBarry Smith PetscCall(MatAssemblyBegin(Prols[iterIdx], MAT_FINAL_ASSEMBLY));
2758be712e4SBarry Smith PetscCall(MatAssemblyEnd(Prols[iterIdx], MAT_FINAL_ASSEMBLY));
2768be712e4SBarry Smith /* project to make new graph for next MIS, skip if last */
2778be712e4SBarry Smith if (iterIdx < misk - 1) {
2788be712e4SBarry Smith Mat new_mat;
279fb842aefSJose E. Roman PetscCall(MatPtAP(cMat, Prols[iterIdx], MAT_INITIAL_MATRIX, PETSC_DETERMINE, &new_mat));
2808be712e4SBarry Smith PetscCall(MatDestroy(&cMat));
2818be712e4SBarry Smith cMat = new_mat; // next iter
2828be712e4SBarry Smith } else if (cMat != Gmat) PetscCall(MatDestroy(&cMat));
2838be712e4SBarry Smith // cleanup
2848be712e4SBarry Smith PetscCall(PetscCDDestroy(agg_lists));
2858be712e4SBarry Smith } /* MIS-k iteration */
2868be712e4SBarry Smith /* make total prolongator Rtot = P_0 * P_1 * ... */
2878be712e4SBarry Smith Rtot = Prols[misk - 1]; // compose P then transpose to get R
2888be712e4SBarry Smith for (PetscInt iterIdx = misk - 1; iterIdx > 0; iterIdx--) {
2898be712e4SBarry Smith Mat P;
290c376f201SBarry Smith
291fb842aefSJose E. Roman PetscCall(MatMatMult(Prols[iterIdx - 1], Rtot, MAT_INITIAL_MATRIX, PETSC_CURRENT, &P));
2928be712e4SBarry Smith PetscCall(MatDestroy(&Prols[iterIdx - 1]));
2938be712e4SBarry Smith PetscCall(MatDestroy(&Rtot));
2948be712e4SBarry Smith Rtot = P;
2958be712e4SBarry Smith }
2968be712e4SBarry Smith PetscCall(MatTranspose(Rtot, MAT_INPLACE_MATRIX, &Rtot)); // R now
2978be712e4SBarry Smith PetscCall(MatViewFromOptions(Rtot, NULL, "-misk_aggregation_view"));
2988be712e4SBarry Smith /* make aggregates with Rtot - could use Rtot directly in theory but have to go through the aggregate list data structure */
2998be712e4SBarry Smith {
3008be712e4SBarry Smith PetscInt Istart, Iend, ncols, NN, MM, jj = 0, max_osz = 0;
3018be712e4SBarry Smith const PetscInt nloc = Gmat->rmap->n;
3028be712e4SBarry Smith PetscCoarsenData *agg_lists;
3038be712e4SBarry Smith Mat mat;
304c376f201SBarry Smith
3058be712e4SBarry Smith PetscCall(PetscCDCreate(nloc, &agg_lists));
3068be712e4SBarry Smith *a_locals_llist = agg_lists; // return
3078be712e4SBarry Smith PetscCall(MatGetOwnershipRange(Rtot, &Istart, &Iend));
308c376f201SBarry Smith for (PetscInt grow = Istart, lid = 0; grow < Iend; grow++, lid++) {
3098be712e4SBarry Smith const PetscInt *idx;
310c376f201SBarry Smith
3118be712e4SBarry Smith PetscCall(MatGetRow(Rtot, grow, &ncols, &idx, NULL));
312c376f201SBarry Smith for (PetscInt jj = 0; jj < ncols; jj++) {
3138be712e4SBarry Smith PetscInt gcol = idx[jj];
314c376f201SBarry Smith
3158be712e4SBarry Smith PetscCall(PetscCDAppendID(agg_lists, lid, gcol)); // local row, global column
3168be712e4SBarry Smith }
3178be712e4SBarry Smith PetscCall(MatRestoreRow(Rtot, grow, &ncols, &idx, NULL));
3188be712e4SBarry Smith }
3198be712e4SBarry Smith PetscCall(MatDestroy(&Rtot));
3208be712e4SBarry Smith
3218be712e4SBarry Smith /* make fake matrix, get largest nnz */
322c376f201SBarry Smith for (PetscInt lid = 0; lid < nloc; lid++) {
3238be712e4SBarry Smith PetscCall(PetscCDCountAt(agg_lists, lid, &jj));
3248be712e4SBarry Smith if (jj > max_osz) max_osz = jj;
3258be712e4SBarry Smith }
3268be712e4SBarry Smith PetscCall(MatGetSize(Gmat, &MM, &NN));
3278be712e4SBarry Smith if (max_osz > MM - nloc) max_osz = MM - nloc;
3288be712e4SBarry Smith PetscCall(MatGetOwnershipRange(Gmat, &Istart, NULL));
3298be712e4SBarry Smith /* matrix of ghost adj for square graph */
3308be712e4SBarry Smith PetscCall(MatCreateAIJ(comm, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE, 0, NULL, max_osz, NULL, &mat));
3318be712e4SBarry Smith for (PetscInt lid = 0, gidi = Istart; lid < nloc; lid++, gidi++) {
3328be712e4SBarry Smith PetscCDIntNd *pos;
333c376f201SBarry Smith
3348be712e4SBarry Smith PetscCall(PetscCDGetHeadPos(agg_lists, lid, &pos));
3358be712e4SBarry Smith while (pos) {
3368be712e4SBarry Smith PetscInt gidj;
337c376f201SBarry Smith
3388be712e4SBarry Smith PetscCall(PetscCDIntNdGetID(pos, &gidj));
3398be712e4SBarry Smith PetscCall(PetscCDGetNextPos(agg_lists, lid, &pos));
3408be712e4SBarry Smith if (gidj < Istart || gidj >= Istart + nloc) PetscCall(MatSetValues(mat, 1, &gidi, 1, &gidj, &one, ADD_VALUES));
3418be712e4SBarry Smith }
3428be712e4SBarry Smith }
3438be712e4SBarry Smith PetscCall(MatAssemblyBegin(mat, MAT_FINAL_ASSEMBLY));
3448be712e4SBarry Smith PetscCall(MatAssemblyEnd(mat, MAT_FINAL_ASSEMBLY));
3458be712e4SBarry Smith PetscCall(PetscCDSetMat(agg_lists, mat));
3468be712e4SBarry Smith }
3478be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3488be712e4SBarry Smith }
3498be712e4SBarry Smith
3508be712e4SBarry Smith /*
3518be712e4SBarry Smith Distance k MIS. k is in 'subctx'
3528be712e4SBarry Smith */
MatCoarsenApply_MISK(MatCoarsen coarse)3538be712e4SBarry Smith static PetscErrorCode MatCoarsenApply_MISK(MatCoarsen coarse)
3548be712e4SBarry Smith {
3558be712e4SBarry Smith Mat mat = coarse->graph;
3568be712e4SBarry Smith PetscInt k;
3578be712e4SBarry Smith
3588be712e4SBarry Smith PetscFunctionBegin;
3598be712e4SBarry Smith PetscCall(MatCoarsenMISKGetDistance(coarse, &k));
360835f2295SStefano Zampini PetscCheck(k > 0, PETSC_COMM_SELF, PETSC_ERR_SUP, "too few levels: %" PetscInt_FMT, k);
3618be712e4SBarry Smith if (!coarse->perm) {
3628be712e4SBarry Smith IS perm;
3638be712e4SBarry Smith PetscInt n, m;
3648be712e4SBarry Smith
3658be712e4SBarry Smith PetscCall(MatGetLocalSize(mat, &m, &n));
3668be712e4SBarry Smith PetscCall(ISCreateStride(PetscObjectComm((PetscObject)mat), m, 0, 1, &perm));
367835f2295SStefano Zampini PetscCall(MatCoarsenApply_MISK_private(perm, k, mat, &coarse->agg_lists));
3688be712e4SBarry Smith PetscCall(ISDestroy(&perm));
3698be712e4SBarry Smith } else {
370835f2295SStefano Zampini PetscCall(MatCoarsenApply_MISK_private(coarse->perm, k, mat, &coarse->agg_lists));
3718be712e4SBarry Smith }
3728be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3738be712e4SBarry Smith }
3748be712e4SBarry Smith
MatCoarsenView_MISK(MatCoarsen coarse,PetscViewer viewer)3758be712e4SBarry Smith static PetscErrorCode MatCoarsenView_MISK(MatCoarsen coarse, PetscViewer viewer)
3768be712e4SBarry Smith {
3778be712e4SBarry Smith PetscMPIInt rank;
378*9f196a02SMartin Diehl PetscBool isascii;
379e0b7e82fSBarry Smith PetscViewerFormat format;
3808be712e4SBarry Smith
3818be712e4SBarry Smith PetscFunctionBegin;
3828be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)coarse), &rank));
383*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
384e0b7e82fSBarry Smith PetscCall(PetscViewerGetFormat(viewer, &format));
385*9f196a02SMartin Diehl if (isascii && format == PETSC_VIEWER_ASCII_INFO_DETAIL && coarse->agg_lists) {
3868be712e4SBarry Smith PetscCall(PetscViewerASCIIPushSynchronized(viewer));
3878be712e4SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] MISK aggregator\n", rank));
3888be712e4SBarry Smith if (!rank) PetscCall(PetscCoarsenDataView_private(coarse->agg_lists, viewer));
3898be712e4SBarry Smith PetscCall(PetscViewerFlush(viewer));
3908be712e4SBarry Smith PetscCall(PetscViewerASCIIPopSynchronized(viewer));
3918be712e4SBarry Smith }
3928be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3938be712e4SBarry Smith }
3948be712e4SBarry Smith
MatCoarsenSetFromOptions_MISK(MatCoarsen coarse,PetscOptionItems PetscOptionsObject)395ce78bad3SBarry Smith static PetscErrorCode MatCoarsenSetFromOptions_MISK(MatCoarsen coarse, PetscOptionItems PetscOptionsObject)
3968be712e4SBarry Smith {
3978be712e4SBarry Smith PetscInt k = 1;
3988be712e4SBarry Smith PetscBool flg;
3994d86920dSPierre Jolivet
4008be712e4SBarry Smith PetscFunctionBegin;
4018be712e4SBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "MatCoarsen-MISk options");
4028be712e4SBarry Smith PetscCall(PetscOptionsInt("-mat_coarsen_misk_distance", "k distance for MIS", "", k, &k, &flg));
4038be712e4SBarry Smith if (flg) coarse->subctx = (void *)(size_t)k;
4048be712e4SBarry Smith PetscOptionsHeadEnd();
4058be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4068be712e4SBarry Smith }
4078be712e4SBarry Smith
4088be712e4SBarry Smith /*MC
4098be712e4SBarry Smith MATCOARSENMISK - A coarsener that uses MISK, a simple greedy coarsener
4108be712e4SBarry Smith
4118be712e4SBarry Smith Level: beginner
4128be712e4SBarry Smith
4138be712e4SBarry Smith Options Database Key:
4148be712e4SBarry Smith . -mat_coarsen_misk_distance <k> - distance for MIS
4158be712e4SBarry Smith
416e0b7e82fSBarry Smith Note:
417e0b7e82fSBarry Smith When the coarsening is used inside `PCGAMG` then the options database key is `-pc_gamg_mat_coarsen_misk_distance`
418e0b7e82fSBarry Smith
419c376f201SBarry Smith .seealso: `MatCoarsen`, `MatCoarsenMISKSetDistance()`, `MatCoarsenApply()`, `MatCoarsenSetType()`, `MatCoarsenType`, `MatCoarsenCreate()`, `MATCOARSENHEM`, `MATCOARSENMIS`
4208be712e4SBarry Smith M*/
4218be712e4SBarry Smith
MatCoarsenCreate_MISK(MatCoarsen coarse)4228be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatCoarsenCreate_MISK(MatCoarsen coarse)
4238be712e4SBarry Smith {
4248be712e4SBarry Smith PetscFunctionBegin;
4258be712e4SBarry Smith coarse->ops->apply = MatCoarsenApply_MISK;
4268be712e4SBarry Smith coarse->ops->view = MatCoarsenView_MISK;
4278be712e4SBarry Smith coarse->subctx = (void *)(size_t)1;
4288be712e4SBarry Smith coarse->ops->setfromoptions = MatCoarsenSetFromOptions_MISK;
4298be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4308be712e4SBarry Smith }
4318be712e4SBarry Smith
4328be712e4SBarry Smith /*@
4338be712e4SBarry Smith MatCoarsenMISKSetDistance - the distance to be used by MISK
4348be712e4SBarry Smith
4358be712e4SBarry Smith Collective
4368be712e4SBarry Smith
4378be712e4SBarry Smith Input Parameters:
4388be712e4SBarry Smith + crs - the coarsen
4398be712e4SBarry Smith - k - the distance
4408be712e4SBarry Smith
4418be712e4SBarry Smith Options Database Key:
4428be712e4SBarry Smith . -mat_coarsen_misk_distance <k> - distance for MIS
4438be712e4SBarry Smith
4448be712e4SBarry Smith Level: advanced
4458be712e4SBarry Smith
446e0b7e82fSBarry Smith Note:
447e0b7e82fSBarry Smith When the coarsening is used inside `PCGAMG` then the options database key is `-pc_gamg_mat_coarsen_misk_distance`
448e0b7e82fSBarry Smith
4490241b274SPierre Jolivet .seealso: `MATCOARSENMISK`, `MatCoarsen`, `MatCoarsenSetFromOptions()`, `MatCoarsenSetType()`, `MatCoarsenRegister()`, `MatCoarsenCreate()`,
4508be712e4SBarry Smith `MatCoarsenDestroy()`, `MatCoarsenSetAdjacency()`, `MatCoarsenMISKGetDistance()`
4518be712e4SBarry Smith `MatCoarsenGetData()`
4528be712e4SBarry Smith @*/
MatCoarsenMISKSetDistance(MatCoarsen crs,PetscInt k)4538be712e4SBarry Smith PetscErrorCode MatCoarsenMISKSetDistance(MatCoarsen crs, PetscInt k)
4548be712e4SBarry Smith {
4558be712e4SBarry Smith PetscFunctionBegin;
4568be712e4SBarry Smith crs->subctx = (void *)(size_t)k;
4578be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4588be712e4SBarry Smith }
4598be712e4SBarry Smith
4608be712e4SBarry Smith /*@
4618be712e4SBarry Smith MatCoarsenMISKGetDistance - gets the distance to be used by MISK
4628be712e4SBarry Smith
4638be712e4SBarry Smith Collective
4648be712e4SBarry Smith
4658be712e4SBarry Smith Input Parameter:
4668be712e4SBarry Smith . crs - the coarsen
4678be712e4SBarry Smith
4688be712e4SBarry Smith Output Parameter:
4698be712e4SBarry Smith . k - the distance
4708be712e4SBarry Smith
4718be712e4SBarry Smith Level: advanced
4728be712e4SBarry Smith
4730241b274SPierre Jolivet .seealso: `MATCOARSENMISK`, `MatCoarsen`, `MatCoarsenSetFromOptions()`, `MatCoarsenSetType()`,
4748be712e4SBarry Smith `MatCoarsenRegister()`, `MatCoarsenCreate()`, `MatCoarsenDestroy()`,
4758be712e4SBarry Smith `MatCoarsenSetAdjacency()`, `MatCoarsenGetData()`
4768be712e4SBarry Smith @*/
MatCoarsenMISKGetDistance(MatCoarsen crs,PetscInt * k)4778be712e4SBarry Smith PetscErrorCode MatCoarsenMISKGetDistance(MatCoarsen crs, PetscInt *k)
4788be712e4SBarry Smith {
4798be712e4SBarry Smith PetscFunctionBegin;
4808be712e4SBarry Smith *k = (PetscInt)(size_t)crs->subctx;
4818be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4828be712e4SBarry Smith }
483