1*8be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/ 2*8be712e4SBarry Smith 3*8be712e4SBarry Smith /* 4*8be712e4SBarry Smith Currently using ParMetis-4.0.2 5*8be712e4SBarry Smith */ 6*8be712e4SBarry Smith 7*8be712e4SBarry Smith #include <parmetis.h> 8*8be712e4SBarry Smith 9*8be712e4SBarry Smith /* 10*8be712e4SBarry Smith The first 5 elements of this structure are the input control array to Metis 11*8be712e4SBarry Smith */ 12*8be712e4SBarry Smith typedef struct { 13*8be712e4SBarry Smith PetscInt cuts; /* number of cuts made (output) */ 14*8be712e4SBarry Smith PetscInt foldfactor; 15*8be712e4SBarry Smith PetscInt parallel; /* use parallel partitioner for coarse problem */ 16*8be712e4SBarry Smith PetscInt indexing; /* 0 indicates C indexing, 1 Fortran */ 17*8be712e4SBarry Smith PetscInt printout; /* indicates if one wishes Metis to print info */ 18*8be712e4SBarry Smith PetscBool repartition; 19*8be712e4SBarry Smith } MatPartitioning_Parmetis; 20*8be712e4SBarry Smith 21*8be712e4SBarry Smith #define PetscCallPARMETIS(n, func) \ 22*8be712e4SBarry Smith do { \ 23*8be712e4SBarry Smith PetscCheck(n != METIS_ERROR_INPUT, PETSC_COMM_SELF, PETSC_ERR_LIB, "ParMETIS error due to wrong inputs and/or options for %s", func); \ 24*8be712e4SBarry Smith PetscCheck(n != METIS_ERROR_MEMORY, PETSC_COMM_SELF, PETSC_ERR_LIB, "ParMETIS error due to insufficient memory in %s", func); \ 25*8be712e4SBarry Smith PetscCheck(n != METIS_ERROR, PETSC_COMM_SELF, PETSC_ERR_LIB, "ParMETIS general error in %s", func); \ 26*8be712e4SBarry Smith } while (0) 27*8be712e4SBarry Smith 28*8be712e4SBarry Smith #define PetscCallParmetis_(name, func, args) \ 29*8be712e4SBarry Smith do { \ 30*8be712e4SBarry Smith PetscStackPushExternal(name); \ 31*8be712e4SBarry Smith int status = func args; \ 32*8be712e4SBarry Smith PetscStackPop; \ 33*8be712e4SBarry Smith PetscCallPARMETIS(status, name); \ 34*8be712e4SBarry Smith } while (0) 35*8be712e4SBarry Smith 36*8be712e4SBarry Smith #define PetscCallParmetis(func, args) PetscCallParmetis_(PetscStringize(func), func, args) 37*8be712e4SBarry Smith 38*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Parmetis_Private(MatPartitioning part, PetscBool useND, PetscBool isImprove, IS *partitioning) 39*8be712e4SBarry Smith { 40*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 41*8be712e4SBarry Smith PetscInt *locals = NULL; 42*8be712e4SBarry Smith Mat mat = part->adj, amat, pmat; 43*8be712e4SBarry Smith PetscBool flg; 44*8be712e4SBarry Smith PetscInt bs = 1; 45*8be712e4SBarry Smith 46*8be712e4SBarry Smith PetscFunctionBegin; 47*8be712e4SBarry Smith PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1); 48*8be712e4SBarry Smith PetscAssertPointer(partitioning, 4); 49*8be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg)); 50*8be712e4SBarry Smith if (flg) { 51*8be712e4SBarry Smith amat = mat; 52*8be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)amat)); 53*8be712e4SBarry Smith } else { 54*8be712e4SBarry Smith /* bs indicates if the converted matrix is "reduced" from the original and hence the 55*8be712e4SBarry Smith resulting partition results need to be stretched to match the original matrix */ 56*8be712e4SBarry Smith PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &amat)); 57*8be712e4SBarry Smith if (amat->rmap->n > 0) bs = mat->rmap->n / amat->rmap->n; 58*8be712e4SBarry Smith } 59*8be712e4SBarry Smith PetscCall(MatMPIAdjCreateNonemptySubcommMat(amat, &pmat)); 60*8be712e4SBarry Smith PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)part))); 61*8be712e4SBarry Smith 62*8be712e4SBarry Smith if (pmat) { 63*8be712e4SBarry Smith MPI_Comm pcomm, comm; 64*8be712e4SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)pmat->data; 65*8be712e4SBarry Smith PetscInt *vtxdist = pmat->rmap->range; 66*8be712e4SBarry Smith PetscInt *xadj = adj->i; 67*8be712e4SBarry Smith PetscInt *adjncy = adj->j; 68*8be712e4SBarry Smith PetscInt *NDorder = NULL; 69*8be712e4SBarry Smith PetscInt itmp = 0, wgtflag = 0, numflag = 0, ncon = part->ncon, nparts = part->n, options[24], i, j; 70*8be712e4SBarry Smith real_t *tpwgts, *ubvec, itr = 0.1; 71*8be712e4SBarry Smith 72*8be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)pmat, &pcomm)); 73*8be712e4SBarry Smith if (PetscDefined(USE_DEBUG)) { 74*8be712e4SBarry Smith /* check that matrix has no diagonal entries */ 75*8be712e4SBarry Smith PetscInt rstart; 76*8be712e4SBarry Smith PetscCall(MatGetOwnershipRange(pmat, &rstart, NULL)); 77*8be712e4SBarry Smith for (i = 0; i < pmat->rmap->n; i++) { 78*8be712e4SBarry Smith for (j = xadj[i]; j < xadj[i + 1]; j++) PetscCheck(adjncy[j] != i + rstart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Row %" PetscInt_FMT " has diagonal entry; Parmetis forbids diagonal entry", i + rstart); 79*8be712e4SBarry Smith } 80*8be712e4SBarry Smith } 81*8be712e4SBarry Smith 82*8be712e4SBarry Smith PetscCall(PetscMalloc1(pmat->rmap->n, &locals)); 83*8be712e4SBarry Smith 84*8be712e4SBarry Smith if (isImprove) { 85*8be712e4SBarry Smith PetscInt i; 86*8be712e4SBarry Smith const PetscInt *part_indices; 87*8be712e4SBarry Smith PetscValidHeaderSpecific(*partitioning, IS_CLASSID, 4); 88*8be712e4SBarry Smith PetscCall(ISGetIndices(*partitioning, &part_indices)); 89*8be712e4SBarry Smith for (i = 0; i < pmat->rmap->n; i++) locals[i] = part_indices[i * bs]; 90*8be712e4SBarry Smith PetscCall(ISRestoreIndices(*partitioning, &part_indices)); 91*8be712e4SBarry Smith PetscCall(ISDestroy(partitioning)); 92*8be712e4SBarry Smith } 93*8be712e4SBarry Smith 94*8be712e4SBarry Smith if (adj->values && part->use_edge_weights && !part->vertex_weights) wgtflag = 1; 95*8be712e4SBarry Smith if (part->vertex_weights && !adj->values) wgtflag = 2; 96*8be712e4SBarry Smith if (part->vertex_weights && adj->values && part->use_edge_weights) wgtflag = 3; 97*8be712e4SBarry Smith 98*8be712e4SBarry Smith if (PetscLogPrintInfo) { 99*8be712e4SBarry Smith itmp = pmetis->printout; 100*8be712e4SBarry Smith pmetis->printout = 127; 101*8be712e4SBarry Smith } 102*8be712e4SBarry Smith PetscCall(PetscMalloc1(ncon * nparts, &tpwgts)); 103*8be712e4SBarry Smith for (i = 0; i < ncon; i++) { 104*8be712e4SBarry Smith for (j = 0; j < nparts; j++) { 105*8be712e4SBarry Smith if (part->part_weights) { 106*8be712e4SBarry Smith tpwgts[i * nparts + j] = part->part_weights[i * nparts + j]; 107*8be712e4SBarry Smith } else { 108*8be712e4SBarry Smith tpwgts[i * nparts + j] = 1. / nparts; 109*8be712e4SBarry Smith } 110*8be712e4SBarry Smith } 111*8be712e4SBarry Smith } 112*8be712e4SBarry Smith PetscCall(PetscMalloc1(ncon, &ubvec)); 113*8be712e4SBarry Smith for (i = 0; i < ncon; i++) ubvec[i] = 1.05; 114*8be712e4SBarry Smith /* This sets the defaults */ 115*8be712e4SBarry Smith options[0] = 1; 116*8be712e4SBarry Smith for (i = 1; i < 24; i++) options[i] = -1; 117*8be712e4SBarry Smith options[1] = 0; /* no verbosity */ 118*8be712e4SBarry Smith options[2] = 0; 119*8be712e4SBarry Smith options[3] = PARMETIS_PSR_COUPLED; /* Seed */ 120*8be712e4SBarry Smith /* Duplicate the communicator to be sure that ParMETIS attribute caching does not interfere with PETSc. */ 121*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_dup(pcomm, &comm)); 122*8be712e4SBarry Smith if (useND) { 123*8be712e4SBarry Smith PetscInt *sizes, *seps, log2size, subd, *level; 124*8be712e4SBarry Smith PetscMPIInt size; 125*8be712e4SBarry Smith idx_t mtype = PARMETIS_MTYPE_GLOBAL, rtype = PARMETIS_SRTYPE_2PHASE, p_nseps = 1, s_nseps = 1; 126*8be712e4SBarry Smith real_t ubfrac = 1.05; 127*8be712e4SBarry Smith 128*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 129*8be712e4SBarry Smith PetscCall(PetscMalloc1(pmat->rmap->n, &NDorder)); 130*8be712e4SBarry Smith PetscCall(PetscMalloc3(2 * size, &sizes, 4 * size, &seps, size, &level)); 131*8be712e4SBarry Smith PetscCallParmetis(ParMETIS_V32_NodeND, ((idx_t *)vtxdist, (idx_t *)xadj, (idx_t *)adjncy, (idx_t *)part->vertex_weights, (idx_t *)&numflag, &mtype, &rtype, &p_nseps, &s_nseps, &ubfrac, NULL /* seed */, NULL /* dbglvl */, (idx_t *)NDorder, (idx_t *)(sizes), &comm)); 132*8be712e4SBarry Smith log2size = PetscLog2Real(size); 133*8be712e4SBarry Smith subd = PetscPowInt(2, log2size); 134*8be712e4SBarry Smith PetscCall(MatPartitioningSizesToSep_Private(subd, sizes, seps, level)); 135*8be712e4SBarry Smith for (i = 0; i < pmat->rmap->n; i++) { 136*8be712e4SBarry Smith PetscInt loc; 137*8be712e4SBarry Smith 138*8be712e4SBarry Smith PetscCall(PetscFindInt(NDorder[i], 2 * subd, seps, &loc)); 139*8be712e4SBarry Smith if (loc < 0) { 140*8be712e4SBarry Smith loc = -(loc + 1); 141*8be712e4SBarry Smith if (loc % 2) { /* part of subdomain */ 142*8be712e4SBarry Smith locals[i] = loc / 2; 143*8be712e4SBarry Smith } else { 144*8be712e4SBarry Smith PetscCall(PetscFindInt(NDorder[i], 2 * (subd - 1), seps + 2 * subd, &loc)); 145*8be712e4SBarry Smith loc = loc < 0 ? -(loc + 1) / 2 : loc / 2; 146*8be712e4SBarry Smith locals[i] = level[loc]; 147*8be712e4SBarry Smith } 148*8be712e4SBarry Smith } else locals[i] = loc / 2; 149*8be712e4SBarry Smith } 150*8be712e4SBarry Smith PetscCall(PetscFree3(sizes, seps, level)); 151*8be712e4SBarry Smith } else { 152*8be712e4SBarry Smith if (pmetis->repartition) { 153*8be712e4SBarry Smith PetscCallParmetis(ParMETIS_V3_AdaptiveRepart, ((idx_t *)vtxdist, (idx_t *)xadj, (idx_t *)adjncy, (idx_t *)part->vertex_weights, (idx_t *)part->vertex_weights, (idx_t *)adj->values, (idx_t *)&wgtflag, (idx_t *)&numflag, (idx_t *)&ncon, (idx_t *)&nparts, tpwgts, ubvec, &itr, (idx_t *)options, 154*8be712e4SBarry Smith (idx_t *)&pmetis->cuts, (idx_t *)locals, &comm)); 155*8be712e4SBarry Smith } else if (isImprove) { 156*8be712e4SBarry Smith PetscCallParmetis(ParMETIS_V3_RefineKway, ((idx_t *)vtxdist, (idx_t *)xadj, (idx_t *)adjncy, (idx_t *)part->vertex_weights, (idx_t *)adj->values, (idx_t *)&wgtflag, (idx_t *)&numflag, (idx_t *)&ncon, (idx_t *)&nparts, tpwgts, ubvec, (idx_t *)options, 157*8be712e4SBarry Smith (idx_t *)&pmetis->cuts, (idx_t *)locals, &comm)); 158*8be712e4SBarry Smith } else { 159*8be712e4SBarry Smith PetscCallParmetis(ParMETIS_V3_PartKway, ((idx_t *)vtxdist, (idx_t *)xadj, (idx_t *)adjncy, (idx_t *)part->vertex_weights, (idx_t *)adj->values, (idx_t *)&wgtflag, (idx_t *)&numflag, (idx_t *)&ncon, (idx_t *)&nparts, tpwgts, ubvec, (idx_t *)options, 160*8be712e4SBarry Smith (idx_t *)&pmetis->cuts, (idx_t *)locals, &comm)); 161*8be712e4SBarry Smith } 162*8be712e4SBarry Smith } 163*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_free(&comm)); 164*8be712e4SBarry Smith 165*8be712e4SBarry Smith PetscCall(PetscFree(tpwgts)); 166*8be712e4SBarry Smith PetscCall(PetscFree(ubvec)); 167*8be712e4SBarry Smith if (PetscLogPrintInfo) pmetis->printout = itmp; 168*8be712e4SBarry Smith 169*8be712e4SBarry Smith if (bs > 1) { 170*8be712e4SBarry Smith PetscInt i, j, *newlocals; 171*8be712e4SBarry Smith PetscCall(PetscMalloc1(bs * pmat->rmap->n, &newlocals)); 172*8be712e4SBarry Smith for (i = 0; i < pmat->rmap->n; i++) { 173*8be712e4SBarry Smith for (j = 0; j < bs; j++) newlocals[bs * i + j] = locals[i]; 174*8be712e4SBarry Smith } 175*8be712e4SBarry Smith PetscCall(PetscFree(locals)); 176*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), bs * pmat->rmap->n, newlocals, PETSC_OWN_POINTER, partitioning)); 177*8be712e4SBarry Smith } else { 178*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), pmat->rmap->n, locals, PETSC_OWN_POINTER, partitioning)); 179*8be712e4SBarry Smith } 180*8be712e4SBarry Smith if (useND) { 181*8be712e4SBarry Smith IS ndis; 182*8be712e4SBarry Smith 183*8be712e4SBarry Smith if (bs > 1) { 184*8be712e4SBarry Smith PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)part), bs, pmat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis)); 185*8be712e4SBarry Smith } else { 186*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), pmat->rmap->n, NDorder, PETSC_OWN_POINTER, &ndis)); 187*8be712e4SBarry Smith } 188*8be712e4SBarry Smith PetscCall(ISSetPermutation(ndis)); 189*8be712e4SBarry Smith PetscCall(PetscObjectCompose((PetscObject)(*partitioning), "_petsc_matpartitioning_ndorder", (PetscObject)ndis)); 190*8be712e4SBarry Smith PetscCall(ISDestroy(&ndis)); 191*8be712e4SBarry Smith } 192*8be712e4SBarry Smith } else { 193*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), 0, NULL, PETSC_COPY_VALUES, partitioning)); 194*8be712e4SBarry Smith if (useND) { 195*8be712e4SBarry Smith IS ndis; 196*8be712e4SBarry Smith 197*8be712e4SBarry Smith if (bs > 1) { 198*8be712e4SBarry Smith PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)part), bs, 0, NULL, PETSC_COPY_VALUES, &ndis)); 199*8be712e4SBarry Smith } else { 200*8be712e4SBarry Smith PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), 0, NULL, PETSC_COPY_VALUES, &ndis)); 201*8be712e4SBarry Smith } 202*8be712e4SBarry Smith PetscCall(ISSetPermutation(ndis)); 203*8be712e4SBarry Smith PetscCall(PetscObjectCompose((PetscObject)(*partitioning), "_petsc_matpartitioning_ndorder", (PetscObject)ndis)); 204*8be712e4SBarry Smith PetscCall(ISDestroy(&ndis)); 205*8be712e4SBarry Smith } 206*8be712e4SBarry Smith } 207*8be712e4SBarry Smith PetscCall(MatDestroy(&pmat)); 208*8be712e4SBarry Smith PetscCall(MatDestroy(&amat)); 209*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 210*8be712e4SBarry Smith } 211*8be712e4SBarry Smith 212*8be712e4SBarry Smith /* 213*8be712e4SBarry Smith Uses the ParMETIS parallel matrix partitioner to compute a nested dissection ordering of the matrix in parallel 214*8be712e4SBarry Smith */ 215*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApplyND_Parmetis(MatPartitioning part, IS *partitioning) 216*8be712e4SBarry Smith { 217*8be712e4SBarry Smith PetscFunctionBegin; 218*8be712e4SBarry Smith PetscCall(MatPartitioningApply_Parmetis_Private(part, PETSC_TRUE, PETSC_FALSE, partitioning)); 219*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 220*8be712e4SBarry Smith } 221*8be712e4SBarry Smith 222*8be712e4SBarry Smith /* 223*8be712e4SBarry Smith Uses the ParMETIS parallel matrix partitioner to partition the matrix in parallel 224*8be712e4SBarry Smith */ 225*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Parmetis(MatPartitioning part, IS *partitioning) 226*8be712e4SBarry Smith { 227*8be712e4SBarry Smith PetscFunctionBegin; 228*8be712e4SBarry Smith PetscCall(MatPartitioningApply_Parmetis_Private(part, PETSC_FALSE, PETSC_FALSE, partitioning)); 229*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 230*8be712e4SBarry Smith } 231*8be712e4SBarry Smith 232*8be712e4SBarry Smith /* 233*8be712e4SBarry Smith Uses the ParMETIS to improve the quality of a partition 234*8be712e4SBarry Smith */ 235*8be712e4SBarry Smith static PetscErrorCode MatPartitioningImprove_Parmetis(MatPartitioning part, IS *partitioning) 236*8be712e4SBarry Smith { 237*8be712e4SBarry Smith PetscFunctionBegin; 238*8be712e4SBarry Smith PetscCall(MatPartitioningApply_Parmetis_Private(part, PETSC_FALSE, PETSC_TRUE, partitioning)); 239*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 240*8be712e4SBarry Smith } 241*8be712e4SBarry Smith 242*8be712e4SBarry Smith static PetscErrorCode MatPartitioningView_Parmetis(MatPartitioning part, PetscViewer viewer) 243*8be712e4SBarry Smith { 244*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 245*8be712e4SBarry Smith PetscMPIInt rank; 246*8be712e4SBarry Smith PetscBool iascii; 247*8be712e4SBarry Smith 248*8be712e4SBarry Smith PetscFunctionBegin; 249*8be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank)); 250*8be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 251*8be712e4SBarry Smith if (iascii) { 252*8be712e4SBarry Smith if (pmetis->parallel == 2) { 253*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Using parallel coarse grid partitioner\n")); 254*8be712e4SBarry Smith } else { 255*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Using sequential coarse grid partitioner\n")); 256*8be712e4SBarry Smith } 257*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Using %" PetscInt_FMT " fold factor\n", pmetis->foldfactor)); 258*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 259*8be712e4SBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d]Number of cuts found %" PetscInt_FMT "\n", rank, pmetis->cuts)); 260*8be712e4SBarry Smith PetscCall(PetscViewerFlush(viewer)); 261*8be712e4SBarry Smith PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 262*8be712e4SBarry Smith } 263*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 264*8be712e4SBarry Smith } 265*8be712e4SBarry Smith 266*8be712e4SBarry Smith /*@ 267*8be712e4SBarry Smith MatPartitioningParmetisSetCoarseSequential - Use the sequential code to 268*8be712e4SBarry Smith do the partitioning of the coarse grid. 269*8be712e4SBarry Smith 270*8be712e4SBarry Smith Logically Collective 271*8be712e4SBarry Smith 272*8be712e4SBarry Smith Input Parameter: 273*8be712e4SBarry Smith . part - the partitioning context 274*8be712e4SBarry Smith 275*8be712e4SBarry Smith Level: advanced 276*8be712e4SBarry Smith 277*8be712e4SBarry Smith .seealso: `MATPARTITIONINGPARMETIS` 278*8be712e4SBarry Smith @*/ 279*8be712e4SBarry Smith PetscErrorCode MatPartitioningParmetisSetCoarseSequential(MatPartitioning part) 280*8be712e4SBarry Smith { 281*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 282*8be712e4SBarry Smith 283*8be712e4SBarry Smith PetscFunctionBegin; 284*8be712e4SBarry Smith pmetis->parallel = 1; 285*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 286*8be712e4SBarry Smith } 287*8be712e4SBarry Smith 288*8be712e4SBarry Smith /*@ 289*8be712e4SBarry Smith MatPartitioningParmetisSetRepartition - Repartition 290*8be712e4SBarry Smith current mesh to rebalance computation. 291*8be712e4SBarry Smith 292*8be712e4SBarry Smith Logically Collective 293*8be712e4SBarry Smith 294*8be712e4SBarry Smith Input Parameter: 295*8be712e4SBarry Smith . part - the partitioning context 296*8be712e4SBarry Smith 297*8be712e4SBarry Smith Level: advanced 298*8be712e4SBarry Smith 299*8be712e4SBarry Smith .seealso: `MATPARTITIONINGPARMETIS` 300*8be712e4SBarry Smith @*/ 301*8be712e4SBarry Smith PetscErrorCode MatPartitioningParmetisSetRepartition(MatPartitioning part) 302*8be712e4SBarry Smith { 303*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 304*8be712e4SBarry Smith 305*8be712e4SBarry Smith PetscFunctionBegin; 306*8be712e4SBarry Smith pmetis->repartition = PETSC_TRUE; 307*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 308*8be712e4SBarry Smith } 309*8be712e4SBarry Smith 310*8be712e4SBarry Smith /*@ 311*8be712e4SBarry Smith MatPartitioningParmetisGetEdgeCut - Returns the number of edge cuts in the vertex partition. 312*8be712e4SBarry Smith 313*8be712e4SBarry Smith Input Parameter: 314*8be712e4SBarry Smith . part - the partitioning context 315*8be712e4SBarry Smith 316*8be712e4SBarry Smith Output Parameter: 317*8be712e4SBarry Smith . cut - the edge cut 318*8be712e4SBarry Smith 319*8be712e4SBarry Smith Level: advanced 320*8be712e4SBarry Smith 321*8be712e4SBarry Smith .seealso: `MATPARTITIONINGPARMETIS` 322*8be712e4SBarry Smith @*/ 323*8be712e4SBarry Smith PetscErrorCode MatPartitioningParmetisGetEdgeCut(MatPartitioning part, PetscInt *cut) 324*8be712e4SBarry Smith { 325*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 326*8be712e4SBarry Smith 327*8be712e4SBarry Smith PetscFunctionBegin; 328*8be712e4SBarry Smith *cut = pmetis->cuts; 329*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 330*8be712e4SBarry Smith } 331*8be712e4SBarry Smith 332*8be712e4SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_Parmetis(MatPartitioning part, PetscOptionItems *PetscOptionsObject) 333*8be712e4SBarry Smith { 334*8be712e4SBarry Smith PetscBool flag = PETSC_FALSE; 335*8be712e4SBarry Smith 336*8be712e4SBarry Smith PetscFunctionBegin; 337*8be712e4SBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Set ParMeTiS partitioning options"); 338*8be712e4SBarry Smith PetscCall(PetscOptionsBool("-mat_partitioning_parmetis_coarse_sequential", "Use sequential coarse partitioner", "MatPartitioningParmetisSetCoarseSequential", flag, &flag, NULL)); 339*8be712e4SBarry Smith if (flag) PetscCall(MatPartitioningParmetisSetCoarseSequential(part)); 340*8be712e4SBarry Smith PetscCall(PetscOptionsBool("-mat_partitioning_parmetis_repartition", "", "MatPartitioningParmetisSetRepartition", flag, &flag, NULL)); 341*8be712e4SBarry Smith if (flag) PetscCall(MatPartitioningParmetisSetRepartition(part)); 342*8be712e4SBarry Smith PetscOptionsHeadEnd(); 343*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 344*8be712e4SBarry Smith } 345*8be712e4SBarry Smith 346*8be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_Parmetis(MatPartitioning part) 347*8be712e4SBarry Smith { 348*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis *)part->data; 349*8be712e4SBarry Smith 350*8be712e4SBarry Smith PetscFunctionBegin; 351*8be712e4SBarry Smith PetscCall(PetscFree(pmetis)); 352*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 353*8be712e4SBarry Smith } 354*8be712e4SBarry Smith 355*8be712e4SBarry Smith /*MC 356*8be712e4SBarry Smith MATPARTITIONINGPARMETIS - Creates a partitioning context via the external package PARMETIS. 357*8be712e4SBarry Smith 358*8be712e4SBarry Smith Collective 359*8be712e4SBarry Smith 360*8be712e4SBarry Smith Input Parameter: 361*8be712e4SBarry Smith . part - the partitioning context 362*8be712e4SBarry Smith 363*8be712e4SBarry Smith Options Database Key: 364*8be712e4SBarry Smith . -mat_partitioning_parmetis_coarse_sequential - use sequential PARMETIS coarse partitioner 365*8be712e4SBarry Smith 366*8be712e4SBarry Smith Level: beginner 367*8be712e4SBarry Smith 368*8be712e4SBarry Smith Note: 369*8be712e4SBarry Smith See https://www-users.cs.umn.edu/~karypis/metis/ 370*8be712e4SBarry Smith 371*8be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningParmetisSetCoarseSequential()`, `MatPartitioningParmetisSetRepartition()`, 372*8be712e4SBarry Smith `MatPartitioningParmetisGetEdgeCut()` 373*8be712e4SBarry Smith M*/ 374*8be712e4SBarry Smith 375*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Parmetis(MatPartitioning part) 376*8be712e4SBarry Smith { 377*8be712e4SBarry Smith MatPartitioning_Parmetis *pmetis; 378*8be712e4SBarry Smith 379*8be712e4SBarry Smith PetscFunctionBegin; 380*8be712e4SBarry Smith PetscCall(PetscNew(&pmetis)); 381*8be712e4SBarry Smith part->data = (void *)pmetis; 382*8be712e4SBarry Smith 383*8be712e4SBarry Smith pmetis->cuts = 0; /* output variable */ 384*8be712e4SBarry Smith pmetis->foldfactor = 150; /*folding factor */ 385*8be712e4SBarry Smith pmetis->parallel = 2; /* use parallel partitioner for coarse grid */ 386*8be712e4SBarry Smith pmetis->indexing = 0; /* index numbering starts from 0 */ 387*8be712e4SBarry Smith pmetis->printout = 0; /* print no output while running */ 388*8be712e4SBarry Smith pmetis->repartition = PETSC_FALSE; 389*8be712e4SBarry Smith 390*8be712e4SBarry Smith part->ops->apply = MatPartitioningApply_Parmetis; 391*8be712e4SBarry Smith part->ops->applynd = MatPartitioningApplyND_Parmetis; 392*8be712e4SBarry Smith part->ops->improve = MatPartitioningImprove_Parmetis; 393*8be712e4SBarry Smith part->ops->view = MatPartitioningView_Parmetis; 394*8be712e4SBarry Smith part->ops->destroy = MatPartitioningDestroy_Parmetis; 395*8be712e4SBarry Smith part->ops->setfromoptions = MatPartitioningSetFromOptions_Parmetis; 396*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 397*8be712e4SBarry Smith } 398*8be712e4SBarry Smith 399*8be712e4SBarry Smith /*@ 400*8be712e4SBarry Smith MatMeshToCellGraph - Convert a mesh to a cell graph. 401*8be712e4SBarry Smith 402*8be712e4SBarry Smith Collective 403*8be712e4SBarry Smith 404*8be712e4SBarry Smith Input Parameters: 405*8be712e4SBarry Smith + mesh - the graph that represents the coupling of the vertices of the mesh 406*8be712e4SBarry Smith - ncommonnodes - mesh elements that share this number of common nodes are considered neighbors, use 2 for triangles and 407*8be712e4SBarry Smith quadrilaterials, 3 for tetrahedrals and 4 for hexahedrals 408*8be712e4SBarry Smith 409*8be712e4SBarry Smith Output Parameter: 410*8be712e4SBarry Smith . dual - the dual graph 411*8be712e4SBarry Smith 412*8be712e4SBarry Smith Level: advanced 413*8be712e4SBarry Smith 414*8be712e4SBarry Smith Notes: 415*8be712e4SBarry Smith Uses the ParMETIS package to convert a `Mat` that represents coupling of vertices of a mesh 416*8be712e4SBarry Smith to a `Mat` the represents the graph of the coupling between cells (the "dual" graph) and is 417*8be712e4SBarry Smith suitable for partitioning with the `MatPartitioning` object. Use this to partition cells of a 418*8be712e4SBarry Smith mesh. 419*8be712e4SBarry Smith 420*8be712e4SBarry Smith Currently requires ParMetis to be installed and uses ParMETIS_V3_Mesh2Dual() 421*8be712e4SBarry Smith 422*8be712e4SBarry Smith Each row of the mesh object represents a single cell in the mesh. For triangles it has 3 entries, quadrilaterials 4 entries, 423*8be712e4SBarry Smith tetrahedrals 4 entries and hexahedrals 8 entries. You can mix triangles and quadrilaterals in the same mesh, but cannot 424*8be712e4SBarry Smith mix tetrahedrals and hexahedrals 425*8be712e4SBarry Smith The columns of each row of the `Mat` mesh are the global vertex numbers of the vertices of that row's cell. 426*8be712e4SBarry Smith The number of rows in mesh is number of cells, the number of columns is the number of vertices. 427*8be712e4SBarry Smith 428*8be712e4SBarry Smith .seealso: `MatCreateMPIAdj()`, `MatPartitioningCreate()` 429*8be712e4SBarry Smith @*/ 430*8be712e4SBarry Smith PetscErrorCode MatMeshToCellGraph(Mat mesh, PetscInt ncommonnodes, Mat *dual) 431*8be712e4SBarry Smith { 432*8be712e4SBarry Smith PetscInt *newxadj, *newadjncy; 433*8be712e4SBarry Smith PetscInt numflag = 0; 434*8be712e4SBarry Smith Mat_MPIAdj *adj = (Mat_MPIAdj *)mesh->data, *newadj; 435*8be712e4SBarry Smith PetscBool flg; 436*8be712e4SBarry Smith MPI_Comm comm; 437*8be712e4SBarry Smith 438*8be712e4SBarry Smith PetscFunctionBegin; 439*8be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)mesh, MATMPIADJ, &flg)); 440*8be712e4SBarry Smith PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "Must use MPIAdj matrix type"); 441*8be712e4SBarry Smith 442*8be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)mesh, &comm)); 443*8be712e4SBarry Smith PetscCallParmetis(ParMETIS_V3_Mesh2Dual, ((idx_t *)mesh->rmap->range, (idx_t *)adj->i, (idx_t *)adj->j, (idx_t *)&numflag, (idx_t *)&ncommonnodes, (idx_t **)&newxadj, (idx_t **)&newadjncy, &comm)); 444*8be712e4SBarry Smith PetscCall(MatCreateMPIAdj(PetscObjectComm((PetscObject)mesh), mesh->rmap->n, mesh->rmap->N, newxadj, newadjncy, NULL, dual)); 445*8be712e4SBarry Smith newadj = (Mat_MPIAdj *)(*dual)->data; 446*8be712e4SBarry Smith 447*8be712e4SBarry Smith newadj->freeaijwithfree = PETSC_TRUE; /* signal the matrix should be freed with system free since space was allocated by ParMETIS */ 448*8be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 449*8be712e4SBarry Smith } 450