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(...) \ 309371c9d4SSatish Balay do { PetscCheck(!(__VA_ARGS__), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error calling PT-Scotch library"); } while (0) 31abe9303eSLisandro Dalcin 329371c9d4SSatish Balay static int PTScotch_Strategy(PetscInt strategy) { 33abe9303eSLisandro Dalcin switch (strategy) { 34abe9303eSLisandro Dalcin case 0: return SCOTCH_STRATDEFAULT; 35abe9303eSLisandro Dalcin case 1: return SCOTCH_STRATQUALITY; 36abe9303eSLisandro Dalcin case 2: return SCOTCH_STRATSPEED; 37abe9303eSLisandro Dalcin case 3: return SCOTCH_STRATBALANCE; 38abe9303eSLisandro Dalcin case 4: return SCOTCH_STRATSAFETY; 39abe9303eSLisandro Dalcin case 5: return SCOTCH_STRATSCALABILITY; 40abe9303eSLisandro Dalcin case 6: return SCOTCH_STRATRECURSIVE; 41abe9303eSLisandro Dalcin case 7: return SCOTCH_STRATREMAP; 42abe9303eSLisandro Dalcin default: return SCOTCH_STRATDEFAULT; 43abe9303eSLisandro Dalcin } 44abe9303eSLisandro Dalcin } 45abe9303eSLisandro Dalcin 469371c9d4SSatish Balay 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[]) { 47d7cc930eSLisandro Dalcin SCOTCH_Arch archdat; 48abe9303eSLisandro Dalcin SCOTCH_Graph grafdat; 49abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 50abe9303eSLisandro Dalcin SCOTCH_Num vertnbr = n; 51abe9303eSLisandro Dalcin SCOTCH_Num edgenbr = xadj[n]; 52abe9303eSLisandro Dalcin SCOTCH_Num *velotab = vtxwgt; 53abe9303eSLisandro Dalcin SCOTCH_Num *edlotab = adjwgt; 54abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 55abe9303eSLisandro Dalcin double kbalval = imbalance; 56abe9303eSLisandro Dalcin 57abe9303eSLisandro Dalcin PetscFunctionBegin; 58abe9303eSLisandro Dalcin { 59abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", NULL, "3.13", "Use -petscpartitioner_use_vertex_weights")); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_ptscotch_vertex_weight", &flg, NULL)); 62abe9303eSLisandro Dalcin if (!flg) velotab = NULL; 63abe9303eSLisandro Dalcin } 649566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphInit(&grafdat)); 659566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphBuild(&grafdat, 0, vertnbr, xadj, xadj + 1, velotab, NULL, edgenbr, adjncy, edlotab)); 669566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat)); 679566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratGraphMapBuild(&stradat, flagval, nparts, kbalval)); 689566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archInit(&archdat)); 69d7cc930eSLisandro Dalcin if (tpart) { 709566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, nparts, tpart)); 71d7cc930eSLisandro Dalcin } else { 729566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, nparts)); 73d7cc930eSLisandro Dalcin } 749566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_graphMap(&grafdat, &archdat, &stradat, part)); 75abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 76abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 77abe9303eSLisandro Dalcin SCOTCH_graphExit(&grafdat); 78abe9303eSLisandro Dalcin PetscFunctionReturn(0); 79abe9303eSLisandro Dalcin } 80abe9303eSLisandro Dalcin 819371c9d4SSatish Balay 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) { 82abe9303eSLisandro Dalcin PetscMPIInt procglbnbr; 83abe9303eSLisandro Dalcin PetscMPIInt proclocnum; 84abe9303eSLisandro Dalcin SCOTCH_Arch archdat; 85abe9303eSLisandro Dalcin SCOTCH_Dgraph grafdat; 86abe9303eSLisandro Dalcin SCOTCH_Dmapping mappdat; 87abe9303eSLisandro Dalcin SCOTCH_Strat stradat; 88abe9303eSLisandro Dalcin SCOTCH_Num vertlocnbr; 89abe9303eSLisandro Dalcin SCOTCH_Num edgelocnbr; 90abe9303eSLisandro Dalcin SCOTCH_Num *veloloctab = vtxwgt; 91abe9303eSLisandro Dalcin SCOTCH_Num *edloloctab = adjwgt; 92abe9303eSLisandro Dalcin SCOTCH_Num flagval = strategy; 93abe9303eSLisandro Dalcin double kbalval = imbalance; 94abe9303eSLisandro Dalcin 95abe9303eSLisandro Dalcin PetscFunctionBegin; 96abe9303eSLisandro Dalcin { 97abe9303eSLisandro Dalcin PetscBool flg = PETSC_TRUE; 989566063dSJacob Faibussowitsch PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", NULL, "3.13", "Use -petscpartitioner_use_vertex_weights")); 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_ptscotch_vertex_weight", &flg, NULL)); 100abe9303eSLisandro Dalcin if (!flg) veloloctab = NULL; 101abe9303eSLisandro Dalcin } 1029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &procglbnbr)); 1039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &proclocnum)); 104abe9303eSLisandro Dalcin vertlocnbr = vtxdist[proclocnum + 1] - vtxdist[proclocnum]; 105abe9303eSLisandro Dalcin edgelocnbr = xadj[vertlocnbr]; 106abe9303eSLisandro Dalcin 1079566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphInit(&grafdat, comm)); 1089566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphBuild(&grafdat, 0, vertlocnbr, vertlocnbr, xadj, xadj + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adjncy, NULL, edloloctab)); 1099566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat)); 1109566063dSJacob Faibussowitsch PetscCall(SCOTCH_stratDgraphMapBuild(&stradat, flagval, procglbnbr, nparts, kbalval)); 1119566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archInit(&archdat)); 112abe9303eSLisandro Dalcin if (tpart) { /* target partition weights */ 1139566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, nparts, tpart)); 114abe9303eSLisandro Dalcin } else { 1159566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, nparts)); 116abe9303eSLisandro Dalcin } 1179566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphMapInit(&grafdat, &mappdat, &archdat, part)); 1189566063dSJacob Faibussowitsch PetscCallPTSCOTCH(SCOTCH_dgraphMapCompute(&grafdat, &mappdat, &stradat)); 119abe9303eSLisandro Dalcin SCOTCH_dgraphMapExit(&grafdat, &mappdat); 120abe9303eSLisandro Dalcin SCOTCH_archExit(&archdat); 121abe9303eSLisandro Dalcin SCOTCH_stratExit(&stradat); 122abe9303eSLisandro Dalcin SCOTCH_dgraphExit(&grafdat); 123abe9303eSLisandro Dalcin PetscFunctionReturn(0); 124abe9303eSLisandro Dalcin } 125abe9303eSLisandro Dalcin 126abe9303eSLisandro Dalcin #endif /* PETSC_HAVE_PTSCOTCH */ 127abe9303eSLisandro Dalcin 1289371c9d4SSatish Balay static const char *const PTScotchStrategyList[] = {"DEFAULT", "QUALITY", "SPEED", "BALANCE", "SAFETY", "SCALABILITY", "RECURSIVE", "REMAP"}; 129abe9303eSLisandro Dalcin 1309371c9d4SSatish Balay static PetscErrorCode PetscPartitionerDestroy_PTScotch(PetscPartitioner part) { 131abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 132abe9303eSLisandro Dalcin 133abe9303eSLisandro Dalcin PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&p->pcomm)); 1359566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 136abe9303eSLisandro Dalcin PetscFunctionReturn(0); 137abe9303eSLisandro Dalcin } 138abe9303eSLisandro Dalcin 1399371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_PTScotch_ASCII(PetscPartitioner part, PetscViewer viewer) { 140abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 141abe9303eSLisandro Dalcin 142abe9303eSLisandro Dalcin PetscFunctionBegin; 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "using partitioning strategy %s\n", PTScotchStrategyList[p->strategy])); 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "using load imbalance ratio %g\n", (double)p->imbalance)); 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 147abe9303eSLisandro Dalcin PetscFunctionReturn(0); 148abe9303eSLisandro Dalcin } 149abe9303eSLisandro Dalcin 1509371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_PTScotch(PetscPartitioner part, PetscViewer viewer) { 151abe9303eSLisandro Dalcin PetscBool iascii; 152abe9303eSLisandro Dalcin 153abe9303eSLisandro Dalcin PetscFunctionBegin; 154abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 155abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1569566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 1579566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscPartitionerView_PTScotch_ASCII(part, viewer)); 158abe9303eSLisandro Dalcin PetscFunctionReturn(0); 159abe9303eSLisandro Dalcin } 160abe9303eSLisandro Dalcin 1619371c9d4SSatish Balay static PetscErrorCode PetscPartitionerSetFromOptions_PTScotch(PetscPartitioner part, PetscOptionItems *PetscOptionsObject) { 162abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data; 163abe9303eSLisandro Dalcin const char *const *slist = PTScotchStrategyList; 164dd39110bSPierre Jolivet PetscInt nlist = PETSC_STATIC_ARRAY_LENGTH(PTScotchStrategyList); 165abe9303eSLisandro Dalcin PetscBool flag; 166abe9303eSLisandro Dalcin 167abe9303eSLisandro Dalcin PetscFunctionBegin; 168d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner PTScotch Options"); 1699566063dSJacob Faibussowitsch PetscCall(PetscOptionsEList("-petscpartitioner_ptscotch_strategy", "Partitioning strategy", "", slist, nlist, slist[p->strategy], &p->strategy, &flag)); 1709566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-petscpartitioner_ptscotch_imbalance", "Load imbalance ratio", "", p->imbalance, &p->imbalance, &flag)); 171d0609cedSBarry Smith PetscOptionsHeadEnd(); 172abe9303eSLisandro Dalcin PetscFunctionReturn(0); 173abe9303eSLisandro Dalcin } 174abe9303eSLisandro Dalcin 1759371c9d4SSatish Balay static PetscErrorCode PetscPartitionerPartition_PTScotch(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) { 176abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PTSCOTCH) 177abe9303eSLisandro Dalcin MPI_Comm comm; 178abe9303eSLisandro Dalcin PetscInt nvtxs = numVertices; /* The number of vertices in full graph */ 179abe9303eSLisandro Dalcin PetscInt *vtxdist; /* Distribution of vertices across processes */ 180abe9303eSLisandro Dalcin PetscInt *xadj = start; /* Start of edge list for each vertex */ 181abe9303eSLisandro Dalcin PetscInt *adjncy = adjacency; /* Edge lists for all vertices */ 182abe9303eSLisandro Dalcin PetscInt *vwgt = NULL; /* Vertex weights */ 183abe9303eSLisandro Dalcin PetscInt *adjwgt = NULL; /* Edge weights */ 184abe9303eSLisandro Dalcin PetscInt v, i, *assignment, *points; 185abe9303eSLisandro Dalcin PetscMPIInt size, rank, p; 186abe9303eSLisandro Dalcin PetscBool hasempty = PETSC_FALSE; 187abe9303eSLisandro Dalcin PetscInt *tpwgts = NULL; 188abe9303eSLisandro Dalcin 189abe9303eSLisandro Dalcin PetscFunctionBegin; 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)part, &comm)); 1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 1929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1939566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(size + 1, &vtxdist, PetscMax(nvtxs, 1), &assignment)); 194abe9303eSLisandro Dalcin /* Calculate vertex distribution */ 195abe9303eSLisandro Dalcin vtxdist[0] = 0; 1969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm)); 197abe9303eSLisandro Dalcin for (p = 2; p <= size; ++p) { 198abe9303eSLisandro Dalcin hasempty = (PetscBool)(hasempty || !vtxdist[p - 1] || !vtxdist[p]); 199abe9303eSLisandro Dalcin vtxdist[p] += vtxdist[p - 1]; 200abe9303eSLisandro Dalcin } 201abe9303eSLisandro Dalcin /* null graph */ 202abe9303eSLisandro Dalcin if (vtxdist[size] == 0) { 2039566063dSJacob Faibussowitsch PetscCall(PetscFree2(vtxdist, assignment)); 2049566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition)); 205abe9303eSLisandro Dalcin PetscFunctionReturn(0); 206abe9303eSLisandro Dalcin } 207abe9303eSLisandro Dalcin 208abe9303eSLisandro Dalcin /* Calculate vertex weights */ 209abe9303eSLisandro Dalcin if (vertSection) { 2109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvtxs, &vwgt)); 2119566063dSJacob Faibussowitsch for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v])); 212abe9303eSLisandro Dalcin } 213abe9303eSLisandro Dalcin 214abe9303eSLisandro Dalcin /* Calculate partition weights */ 215abe9303eSLisandro Dalcin if (targetSection) { 216abe9303eSLisandro Dalcin PetscInt sumw; 217abe9303eSLisandro Dalcin 2189566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &tpwgts)); 219abe9303eSLisandro Dalcin for (p = 0, sumw = 0; p < nparts; ++p) { 2209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection, p, &tpwgts[p])); 221abe9303eSLisandro Dalcin sumw += tpwgts[p]; 222abe9303eSLisandro Dalcin } 2239566063dSJacob Faibussowitsch if (!sumw) PetscCall(PetscFree(tpwgts)); 224abe9303eSLisandro Dalcin } 225abe9303eSLisandro Dalcin 226abe9303eSLisandro Dalcin { 227abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *pts = (PetscPartitioner_PTScotch *)part->data; 228abe9303eSLisandro Dalcin int strat = PTScotch_Strategy(pts->strategy); 229abe9303eSLisandro Dalcin double imbal = (double)pts->imbalance; 230abe9303eSLisandro Dalcin 2319371c9d4SSatish Balay for (p = 0; !vtxdist[p + 1] && p < size; ++p) 2329371c9d4SSatish Balay ; 233abe9303eSLisandro Dalcin if (vtxdist[p + 1] == vtxdist[size]) { 2349566063dSJacob Faibussowitsch if (rank == p) PetscCall(PTScotch_PartGraph_Seq(strat, imbal, nvtxs, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment)); 235abe9303eSLisandro Dalcin } else { 236abe9303eSLisandro Dalcin MPI_Comm pcomm = pts->pcomm; 237abe9303eSLisandro Dalcin 238abe9303eSLisandro Dalcin if (hasempty) { 239abe9303eSLisandro Dalcin PetscInt cnt; 240abe9303eSLisandro Dalcin 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(pts->pcomm, !!nvtxs, rank, &pcomm)); 242abe9303eSLisandro Dalcin for (p = 0, cnt = 0; p < size; p++) { 243abe9303eSLisandro Dalcin if (vtxdist[p + 1] != vtxdist[p]) { 244abe9303eSLisandro Dalcin vtxdist[cnt + 1] = vtxdist[p + 1]; 245abe9303eSLisandro Dalcin cnt++; 246abe9303eSLisandro Dalcin } 247abe9303eSLisandro Dalcin } 248abe9303eSLisandro Dalcin }; 2499566063dSJacob Faibussowitsch if (nvtxs) PetscCall(PTScotch_PartGraph_MPI(strat, imbal, vtxdist, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment, pcomm)); 2509566063dSJacob Faibussowitsch if (hasempty) PetscCallMPI(MPI_Comm_free(&pcomm)); 251abe9303eSLisandro Dalcin } 252abe9303eSLisandro Dalcin } 2539566063dSJacob Faibussowitsch PetscCall(PetscFree(vwgt)); 2549566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts)); 255abe9303eSLisandro Dalcin 256abe9303eSLisandro Dalcin /* Convert to PetscSection+IS */ 2579566063dSJacob Faibussowitsch for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nvtxs, &points)); 259abe9303eSLisandro Dalcin for (p = 0, i = 0; p < nparts; ++p) { 260abe9303eSLisandro Dalcin for (v = 0; v < nvtxs; ++v) { 261abe9303eSLisandro Dalcin if (assignment[v] == p) points[i++] = v; 262abe9303eSLisandro Dalcin } 263abe9303eSLisandro Dalcin } 26463a3b9bcSJacob Faibussowitsch PetscCheck(i == nvtxs, comm, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs); 2659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition)); 266abe9303eSLisandro Dalcin 2679566063dSJacob Faibussowitsch PetscCall(PetscFree2(vtxdist, assignment)); 268abe9303eSLisandro Dalcin PetscFunctionReturn(0); 269abe9303eSLisandro Dalcin #else 270abe9303eSLisandro Dalcin SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-ptscotch."); 271abe9303eSLisandro Dalcin #endif 272abe9303eSLisandro Dalcin } 273abe9303eSLisandro Dalcin 2749371c9d4SSatish Balay static PetscErrorCode PetscPartitionerInitialize_PTScotch(PetscPartitioner part) { 275abe9303eSLisandro Dalcin PetscFunctionBegin; 276abe9303eSLisandro Dalcin part->noGraph = PETSC_FALSE; 277abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_PTScotch; 278abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_PTScotch; 279abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_PTScotch; 280abe9303eSLisandro Dalcin part->ops->setfromoptions = PetscPartitionerSetFromOptions_PTScotch; 281abe9303eSLisandro Dalcin PetscFunctionReturn(0); 282abe9303eSLisandro Dalcin } 283abe9303eSLisandro Dalcin 284abe9303eSLisandro Dalcin /*MC 285abe9303eSLisandro Dalcin PETSCPARTITIONERPTSCOTCH = "ptscotch" - A PetscPartitioner object using the PT-Scotch library 286abe9303eSLisandro Dalcin 287abe9303eSLisandro Dalcin Level: intermediate 288abe9303eSLisandro Dalcin 289abe9303eSLisandro Dalcin Options Database Keys: 290abe9303eSLisandro Dalcin + -petscpartitioner_ptscotch_strategy <string> - PT-Scotch strategy. Choose one of default quality speed balance safety scalability recursive remap 291abe9303eSLisandro Dalcin - -petscpartitioner_ptscotch_imbalance <val> - Load imbalance ratio 292abe9303eSLisandro Dalcin 293abe9303eSLisandro Dalcin Notes: when the graph is on a single process, this partitioner actually uses Scotch and not PT-Scotch 294abe9303eSLisandro Dalcin 295db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()` 296abe9303eSLisandro Dalcin M*/ 297abe9303eSLisandro Dalcin 2989371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_PTScotch(PetscPartitioner part) { 299abe9303eSLisandro Dalcin PetscPartitioner_PTScotch *p; 300abe9303eSLisandro Dalcin 301abe9303eSLisandro Dalcin PetscFunctionBegin; 302abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 303*4dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&p)); 304abe9303eSLisandro Dalcin part->data = p; 305abe9303eSLisandro Dalcin 3069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)part), &p->pcomm)); 307abe9303eSLisandro Dalcin p->strategy = 0; 308abe9303eSLisandro Dalcin p->imbalance = 0.01; 309abe9303eSLisandro Dalcin 3109566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_PTScotch(part)); 3119566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(PTScotchPartitionerCitation, &PTScotchPartitionerCite)); 312abe9303eSLisandro Dalcin PetscFunctionReturn(0); 313abe9303eSLisandro Dalcin } 314