1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/ 2abe9303eSLisandro Dalcin 3abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 4abe9303eSLisandro Dalcin EXTERN_C_BEGIN 5abe9303eSLisandro Dalcin #include <ptscotch.h> 6abe9303eSLisandro Dalcin EXTERN_C_END 7abe9303eSLisandro Dalcin #endif 8abe9303eSLisandro Dalcin 9abe9303eSLisandro Dalcin PetscBool PTScotchPartitionerCite = PETSC_FALSE; 109371c9d4SSatish Balay const char PTScotchPartitionerCitation[] = "@article{PTSCOTCH,\n" 11abe9303eSLisandro Dalcin " author = {C. Chevalier and F. Pellegrini},\n" 12abe9303eSLisandro Dalcin " title = {{PT-SCOTCH}: a tool for efficient parallel graph ordering},\n" 13abe9303eSLisandro Dalcin " journal = {Parallel Computing},\n" 14abe9303eSLisandro Dalcin " volume = {34},\n" 15abe9303eSLisandro Dalcin " number = {6},\n" 16abe9303eSLisandro Dalcin " pages = {318--331},\n" 17abe9303eSLisandro Dalcin " year = {2008},\n" 18abe9303eSLisandro Dalcin " doi = {https://doi.org/10.1016/j.parco.2007.12.001}\n" 19abe9303eSLisandro Dalcin "}\n"; 20abe9303eSLisandro Dalcin 21abe9303eSLisandro Dalcin typedef struct { 22abe9303eSLisandro Dalcin MPI_Comm pcomm; 23abe9303eSLisandro Dalcin PetscInt strategy; 24abe9303eSLisandro Dalcin PetscReal imbalance; 25abe9303eSLisandro Dalcin } PetscPartitioner_PTScotch; 26abe9303eSLisandro Dalcin 27abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 28abe9303eSLisandro Dalcin 299371c9d4SSatish Balay #define PetscCallPTSCOTCH(...) \ 30d71ae5a4SJacob Faibussowitsch do { \ 31d71ae5a4SJacob Faibussowitsch PetscCheck(!(__VA_ARGS__), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error calling PT-Scotch library"); \ 32d71ae5a4SJacob Faibussowitsch } while (0) 33abe9303eSLisandro Dalcin 34d71ae5a4SJacob Faibussowitsch static int PTScotch_Strategy(PetscInt strategy) 35d71ae5a4SJacob Faibussowitsch { 36abe9303eSLisandro Dalcin switch (strategy) { 37d71ae5a4SJacob Faibussowitsch case 0: 38d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATDEFAULT; 39d71ae5a4SJacob Faibussowitsch case 1: 40d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATQUALITY; 41d71ae5a4SJacob Faibussowitsch case 2: 42d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATSPEED; 43d71ae5a4SJacob Faibussowitsch case 3: 44d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATBALANCE; 45d71ae5a4SJacob Faibussowitsch case 4: 46d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATSAFETY; 47d71ae5a4SJacob Faibussowitsch case 5: 48d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATSCALABILITY; 49d71ae5a4SJacob Faibussowitsch case 6: 50d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATRECURSIVE; 51d71ae5a4SJacob Faibussowitsch case 7: 52d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATREMAP; 53d71ae5a4SJacob Faibussowitsch default: 54d71ae5a4SJacob Faibussowitsch return SCOTCH_STRATDEFAULT; 55abe9303eSLisandro Dalcin } 56abe9303eSLisandro Dalcin } 57abe9303eSLisandro Dalcin 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PTScotch_PartGraph_Seq(SCOTCH_Num strategy, double imbalance, SCOTCH_Num n, SCOTCH_Num xadj[], SCOTCH_Num adjncy[], SCOTCH_Num vtxwgt[], SCOTCH_Num adjwgt[], SCOTCH_Num nparts, SCOTCH_Num tpart[], SCOTCH_Num part[]) 59d71ae5a4SJacob Faibussowitsch { 60d7cc930eSLisandro Dalcin SCOTCH_Arch archdat; 61abe9303eSLisandro Dalcin SCOTCH_Graph grafdat; 62abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 63abe9303eSLisandro Dalcin SCOTCH_Num vertnbr = n; 64abe9303eSLisandro Dalcin SCOTCH_Num edgenbr = xadj[n]; 65abe9303eSLisandro Dalcin SCOTCH_Num *velotab = vtxwgt; 66abe9303eSLisandro Dalcin SCOTCH_Num *edlotab = adjwgt; 67abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 68abe9303eSLisandro Dalcin double kbalval = imbalance; 69abe9303eSLisandro Dalcin 70abe9303eSLisandro Dalcin PetscFunctionBegin; 71fd599e24SMatthew G. Knepley if (!n) PetscFunctionReturn(PETSC_SUCCESS); 72abe9303eSLisandro Dalcin { 73abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 744ead3382SBarry Smith PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", "-petscpartitioner_use_vertex_weights", "3.13", NULL)); 754ead3382SBarry Smith /* 764ead3382SBarry Smith Cannot remove the PetscOptionsGetBool() below since the PetscOptionsDeprecatedNoObject() above is called after the non-deprecated version 774ead3382SBarry Smith has already been checked in PetscPartitionerSetFromOptions(). 784ead3382SBarry Smith */ 794ead3382SBarry Smith PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_use_vertex_weight", &flg, NULL)); 80abe9303eSLisandro Dalcin if (!flg) velotab = NULL; 81abe9303eSLisandro Dalcin } 829566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphInit(&grafdat)); 839566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphBuild(&grafdat, 0, vertnbr, xadj, xadj + 1, velotab, NULL, edgenbr, adjncy, edlotab)); 849566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat)); 859566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratGraphMapBuild(&stradat, flagval, nparts, kbalval)); 869566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archInit(&archdat)); 87d7cc930eSLisandro Dalcin if (tpart) { 88fd599e24SMatthew G. Knepley PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, PetscMin(nparts, n), tpart)); 89d7cc930eSLisandro Dalcin } else { 90fd599e24SMatthew G. Knepley PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, PetscMin(nparts, n))); 91d7cc930eSLisandro Dalcin } 929566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphMap(&grafdat, &archdat, &stradat, part)); 93abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 94abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 95abe9303eSLisandro Dalcin SCOTCH_graphExit(&grafdat); 963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97abe9303eSLisandro Dalcin } 98abe9303eSLisandro Dalcin 99d71ae5a4SJacob Faibussowitsch static PetscErrorCode PTScotch_PartGraph_MPI(SCOTCH_Num strategy, double imbalance, SCOTCH_Num vtxdist[], SCOTCH_Num xadj[], SCOTCH_Num adjncy[], SCOTCH_Num vtxwgt[], SCOTCH_Num adjwgt[], SCOTCH_Num nparts, SCOTCH_Num tpart[], SCOTCH_Num part[], MPI_Comm comm) 100d71ae5a4SJacob Faibussowitsch { 101abe9303eSLisandro Dalcin PetscMPIInt procglbnbr; 102abe9303eSLisandro Dalcin PetscMPIInt proclocnum; 103abe9303eSLisandro Dalcin SCOTCH_Arch archdat; 104abe9303eSLisandro Dalcin SCOTCH_Dgraph grafdat; 105abe9303eSLisandro Dalcin SCOTCH_Dmapping mappdat; 106abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 107abe9303eSLisandro Dalcin SCOTCH_Num vertlocnbr; 108abe9303eSLisandro Dalcin SCOTCH_Num edgelocnbr; 109abe9303eSLisandro Dalcin SCOTCH_Num *veloloctab = vtxwgt; 110abe9303eSLisandro Dalcin SCOTCH_Num *edloloctab = adjwgt; 111abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 112abe9303eSLisandro Dalcin double kbalval = imbalance; 113abe9303eSLisandro Dalcin 114abe9303eSLisandro Dalcin PetscFunctionBegin; 115abe9303eSLisandro Dalcin { 116abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 1174ead3382SBarry Smith PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", "-petscpartitioner_use_vertex_weights", "3.13", NULL)); 1184ead3382SBarry Smith /* 1194ead3382SBarry Smith Cannot remove the PetscOptionsGetBool() below since the PetscOptionsDeprecatedNoObject() above is called after the non-deprecated version 1204ead3382SBarry Smith has already been checked in PetscPartitionerSetFromOptions(). 1214ead3382SBarry Smith */ 1224ead3382SBarry Smith PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_use_vertex_weight", &flg, NULL)); 123abe9303eSLisandro Dalcin if (!flg) veloloctab = NULL; 124abe9303eSLisandro Dalcin } 1259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &procglbnbr)); 1269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &proclocnum)); 127abe9303eSLisandro Dalcin vertlocnbr = vtxdist[proclocnum + 1] - vtxdist[proclocnum]; 128abe9303eSLisandro Dalcin edgelocnbr = xadj[vertlocnbr]; 129abe9303eSLisandro Dalcin 1309566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphInit(&grafdat, comm)); 1319566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphBuild(&grafdat, 0, vertlocnbr, vertlocnbr, xadj, xadj + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adjncy, NULL, edloloctab)); 1329566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat)); 1333ba16761SJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratDgraphMapBuild(&stradat, flagval, procglbnbr, nparts, kbalval)); 1349566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archInit(&archdat)); 135abe9303eSLisandro Dalcin if (tpart) { /* target partition weights */ 1369566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, nparts, tpart)); 137abe9303eSLisandro Dalcin } else { 1389566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, nparts)); 139abe9303eSLisandro Dalcin } 1409566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphMapInit(&grafdat, &mappdat, &archdat, part)); 1419566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphMapCompute(&grafdat, &mappdat, &stradat)); 142abe9303eSLisandro Dalcin SCOTCH_dgraphMapExit(&grafdat, &mappdat); 143abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 144abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 145abe9303eSLisandro Dalcin SCOTCH_dgraphExit(&grafdat); 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147abe9303eSLisandro Dalcin } 148abe9303eSLisandro Dalcin 149abe9303eSLisandro Dalcin #endif /* PETSC_HAVE_PTSCOTCH */ 150abe9303eSLisandro Dalcin 1519371c9d4SSatish Balay static const char *const PTScotchStrategyList[] = {"DEFAULT", "QUALITY", "SPEED", "BALANCE", "SAFETY", "SCALABILITY", "RECURSIVE", "REMAP"}; 152abe9303eSLisandro Dalcin 153d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_PTScotch(PetscPartitioner part) 154d71ae5a4SJacob Faibussowitsch { 155abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 156abe9303eSLisandro Dalcin 157abe9303eSLisandro Dalcin PetscFunctionBegin; 1589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&p->pcomm)); 1599566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 1603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161abe9303eSLisandro Dalcin } 162abe9303eSLisandro Dalcin 163d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_PTScotch_ASCII(PetscPartitioner part, PetscViewer viewer) 164d71ae5a4SJacob Faibussowitsch { 165abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 166abe9303eSLisandro Dalcin 167abe9303eSLisandro Dalcin PetscFunctionBegin; 1689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "using partitioning strategy %s\n", PTScotchStrategyList[p->strategy])); 1709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "using load imbalance ratio %g\n", (double)p->imbalance)); 1719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 173abe9303eSLisandro Dalcin } 174abe9303eSLisandro Dalcin 175d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_PTScotch(PetscPartitioner part, PetscViewer viewer) 176d71ae5a4SJacob Faibussowitsch { 177*9f196a02SMartin Diehl PetscBool isascii; 178abe9303eSLisandro Dalcin 179abe9303eSLisandro Dalcin PetscFunctionBegin; 180abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 181abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 182*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 183*9f196a02SMartin Diehl if (isascii) PetscCall(PetscPartitionerView_PTScotch_ASCII(part, viewer)); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185abe9303eSLisandro Dalcin } 186abe9303eSLisandro Dalcin 187ce78bad3SBarry Smith static PetscErrorCode PetscPartitionerSetFromOptions_PTScotch(PetscPartitioner part, PetscOptionItems PetscOptionsObject) 188d71ae5a4SJacob Faibussowitsch { 189abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 190abe9303eSLisandro Dalcin const char *const *slist = PTScotchStrategyList; 191dd39110bSPierre Jolivet PetscInt nlist = PETSC_STATIC_ARRAY_LENGTH(PTScotchStrategyList); 192abe9303eSLisandro Dalcin PetscBool flag; 193abe9303eSLisandro Dalcin 194abe9303eSLisandro Dalcin PetscFunctionBegin; 195d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner PTScotch Options"); 1969566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-petscpartitioner_ptscotch_strategy", "Partitioning strategy", "", slist, nlist, slist[p->strategy], &p->strategy, &flag)); 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-petscpartitioner_ptscotch_imbalance", "Load imbalance ratio", "", p->imbalance, &p->imbalance, &flag)); 198d0609cedSBarry Smith PetscOptionsHeadEnd(); 1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 200abe9303eSLisandro Dalcin } 201abe9303eSLisandro Dalcin 20221c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_PTScotch(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition) 203d71ae5a4SJacob Faibussowitsch { 204abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 205abe9303eSLisandro Dalcin MPI_Comm comm; 206abe9303eSLisandro Dalcin PetscInt nvtxs = numVertices; /* The number of vertices in full graph */ 207abe9303eSLisandro Dalcin PetscInt *vtxdist; /* Distribution of vertices across processes */ 208abe9303eSLisandro Dalcin PetscInt *xadj = start; /* Start of edge list for each vertex */ 209abe9303eSLisandro Dalcin PetscInt *adjncy = adjacency; /* Edge lists for all vertices */ 210abe9303eSLisandro Dalcin PetscInt *vwgt = NULL; /* Vertex weights */ 211abe9303eSLisandro Dalcin PetscInt *adjwgt = NULL; /* Edge weights */ 212abe9303eSLisandro Dalcin PetscInt v, i, *assignment, *points; 213abe9303eSLisandro Dalcin PetscMPIInt size, rank, p; 214abe9303eSLisandro Dalcin PetscBool hasempty = PETSC_FALSE; 215abe9303eSLisandro Dalcin PetscInt *tpwgts = NULL; 216abe9303eSLisandro Dalcin 217abe9303eSLisandro Dalcin PetscFunctionBegin; 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)part, &comm)); 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 2209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 2219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size + 1, &vtxdist, PetscMax(nvtxs, 1), &assignment)); 222abe9303eSLisandro Dalcin /* Calculate vertex distribution */ 223abe9303eSLisandro Dalcin vtxdist[0] = 0; 2249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm)); 225abe9303eSLisandro Dalcin for (p = 2; p <= size; ++p) { 226abe9303eSLisandro Dalcin hasempty = (PetscBool)(hasempty || !vtxdist[p - 1] || !vtxdist[p]); 227abe9303eSLisandro Dalcin vtxdist[p] += vtxdist[p - 1]; 228abe9303eSLisandro Dalcin } 229abe9303eSLisandro Dalcin /* null graph */ 230abe9303eSLisandro Dalcin if (vtxdist[size] == 0) { 2319566063dSJacob Faibussowitsch PetscCall(PetscFree2(vtxdist, assignment)); 2329566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition)); 2333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 234abe9303eSLisandro Dalcin } 235abe9303eSLisandro Dalcin 236abe9303eSLisandro Dalcin /* Calculate vertex weights */ 237abe9303eSLisandro Dalcin if (vertSection) { 2389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvtxs, &vwgt)); 2399566063dSJacob Faibussowitsch for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v])); 240abe9303eSLisandro Dalcin } 24121c92275SMatthew G. Knepley // Weight edges 24221c92275SMatthew G. Knepley if (edgeSection) { 24321c92275SMatthew G. Knepley PetscCall(PetscMalloc1(xadj[nvtxs], &adjwgt)); 24421c92275SMatthew G. Knepley for (PetscInt e = 0; e < xadj[nvtxs]; ++e) PetscCall(PetscSectionGetDof(edgeSection, e, &adjwgt[e])); 24521c92275SMatthew G. Knepley } 246abe9303eSLisandro Dalcin 247abe9303eSLisandro Dalcin /* Calculate partition weights */ 248abe9303eSLisandro Dalcin if (targetSection) { 249abe9303eSLisandro Dalcin PetscInt sumw; 250abe9303eSLisandro Dalcin 2519566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &tpwgts)); 252abe9303eSLisandro Dalcin for (p = 0, sumw = 0; p < nparts; ++p) { 2539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection, p, &tpwgts[p])); 254abe9303eSLisandro Dalcin sumw += tpwgts[p]; 255abe9303eSLisandro Dalcin } 2569566063dSJacob Faibussowitsch if (!sumw) PetscCall(PetscFree(tpwgts)); 257abe9303eSLisandro Dalcin } 258abe9303eSLisandro Dalcin 259abe9303eSLisandro Dalcin { 260abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *pts = (PetscPartitioner_PTScotch *)part->data; 261abe9303eSLisandro Dalcin int strat = PTScotch_Strategy(pts->strategy); 262abe9303eSLisandro Dalcin double imbal = (double)pts->imbalance; 263abe9303eSLisandro Dalcin 264fbccb6d4SPierre Jolivet for (p = 0; !vtxdist[p + 1] && p < size; ++p); 265abe9303eSLisandro Dalcin if (vtxdist[p + 1] == vtxdist[size]) { 2669566063dSJacob Faibussowitsch if (rank == p) PetscCall(PTScotch_PartGraph_Seq(strat, imbal, nvtxs, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment)); 267abe9303eSLisandro Dalcin } else { 268abe9303eSLisandro Dalcin MPI_Comm pcomm = pts->pcomm; 269abe9303eSLisandro Dalcin 270abe9303eSLisandro Dalcin if (hasempty) { 271abe9303eSLisandro Dalcin PetscInt cnt; 272abe9303eSLisandro Dalcin 2739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(pts->pcomm, !!nvtxs, rank, &pcomm)); 274abe9303eSLisandro Dalcin for (p = 0, cnt = 0; p < size; p++) { 275abe9303eSLisandro Dalcin if (vtxdist[p + 1] != vtxdist[p]) { 276abe9303eSLisandro Dalcin vtxdist[cnt + 1] = vtxdist[p + 1]; 277abe9303eSLisandro Dalcin cnt++; 278abe9303eSLisandro Dalcin } 279abe9303eSLisandro Dalcin } 280a8f51744SPierre Jolivet } 2819566063dSJacob Faibussowitsch if (nvtxs) PetscCall(PTScotch_PartGraph_MPI(strat, imbal, vtxdist, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment, pcomm)); 2829566063dSJacob Faibussowitsch if (hasempty) PetscCallMPI(MPI_Comm_free(&pcomm)); 283abe9303eSLisandro Dalcin } 284abe9303eSLisandro Dalcin } 2859566063dSJacob Faibussowitsch PetscCall(PetscFree(vwgt)); 28621c92275SMatthew G. Knepley PetscCall(PetscFree(adjwgt)); 2879566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts)); 288abe9303eSLisandro Dalcin 289abe9303eSLisandro Dalcin /* Convert to PetscSection+IS */ 2909566063dSJacob Faibussowitsch for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1)); 2919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvtxs, &points)); 292abe9303eSLisandro Dalcin for (p = 0, i = 0; p < nparts; ++p) { 293abe9303eSLisandro Dalcin for (v = 0; v < nvtxs; ++v) { 294abe9303eSLisandro Dalcin if (assignment[v] == p) points[i++] = v; 295abe9303eSLisandro Dalcin } 296abe9303eSLisandro Dalcin } 29763a3b9bcSJacob Faibussowitsch PetscCheck(i == nvtxs, comm, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs); 2989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition)); 299abe9303eSLisandro Dalcin 3009566063dSJacob Faibussowitsch PetscCall(PetscFree2(vtxdist, assignment)); 3013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 302abe9303eSLisandro Dalcin #else 303abe9303eSLisandro Dalcin SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-ptscotch."); 304abe9303eSLisandro Dalcin #endif 305abe9303eSLisandro Dalcin } 306abe9303eSLisandro Dalcin 307d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_PTScotch(PetscPartitioner part) 308d71ae5a4SJacob Faibussowitsch { 309abe9303eSLisandro Dalcin PetscFunctionBegin; 310abe9303eSLisandro Dalcin part->noGraph = PETSC_FALSE; 311abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_PTScotch; 312abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_PTScotch; 313abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_PTScotch; 314abe9303eSLisandro Dalcin part->ops->setfromoptions = PetscPartitionerSetFromOptions_PTScotch; 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316abe9303eSLisandro Dalcin } 317abe9303eSLisandro Dalcin 318abe9303eSLisandro Dalcin /*MC 319abe9303eSLisandro Dalcin PETSCPARTITIONERPTSCOTCH = "ptscotch" - A PetscPartitioner object using the PT-Scotch library 320abe9303eSLisandro Dalcin 321abe9303eSLisandro Dalcin Level: intermediate 322abe9303eSLisandro Dalcin 323abe9303eSLisandro Dalcin Options Database Keys: 324abe9303eSLisandro Dalcin + -petscpartitioner_ptscotch_strategy <string> - PT-Scotch strategy. Choose one of default quality speed balance safety scalability recursive remap 325abe9303eSLisandro Dalcin - -petscpartitioner_ptscotch_imbalance <val> - Load imbalance ratio 326abe9303eSLisandro Dalcin 327abe9303eSLisandro Dalcin Notes: when the graph is on a single process, this partitioner actually uses Scotch and not PT-Scotch 328abe9303eSLisandro Dalcin 329db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()` 330abe9303eSLisandro Dalcin M*/ 331abe9303eSLisandro Dalcin 332d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_PTScotch(PetscPartitioner part) 333d71ae5a4SJacob Faibussowitsch { 334abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p; 335abe9303eSLisandro Dalcin 336abe9303eSLisandro Dalcin PetscFunctionBegin; 337abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 3384dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&p)); 339abe9303eSLisandro Dalcin part->data = p; 340abe9303eSLisandro Dalcin 3419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)part), &p->pcomm)); 342abe9303eSLisandro Dalcin p->strategy = 0; 343abe9303eSLisandro Dalcin p->imbalance = 0.01; 344abe9303eSLisandro Dalcin 3459566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_PTScotch(part)); 3469566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(PTScotchPartitionerCitation, &PTScotchPartitionerCite)); 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348abe9303eSLisandro Dalcin } 349