xref: /petsc/src/dm/partitioner/impls/ptscotch/partptscotch.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1 #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
2 
3 #if defined(PETSC_HAVE_PTSCOTCH)
4 EXTERN_C_BEGIN
5   #include <ptscotch.h>
6 EXTERN_C_END
7 #endif
8 
9 PetscBool  PTScotchPartitionerCite       = PETSC_FALSE;
10 const char PTScotchPartitionerCitation[] = "@article{PTSCOTCH,\n"
11                                            "  author  = {C. Chevalier and F. Pellegrini},\n"
12                                            "  title   = {{PT-SCOTCH}: a tool for efficient parallel graph ordering},\n"
13                                            "  journal = {Parallel Computing},\n"
14                                            "  volume  = {34},\n"
15                                            "  number  = {6},\n"
16                                            "  pages   = {318--331},\n"
17                                            "  year    = {2008},\n"
18                                            "  doi     = {https://doi.org/10.1016/j.parco.2007.12.001}\n"
19                                            "}\n";
20 
21 typedef struct {
22   MPI_Comm  pcomm;
23   PetscInt  strategy;
24   PetscReal imbalance;
25 } PetscPartitioner_PTScotch;
26 
27 #if defined(PETSC_HAVE_PTSCOTCH)
28 
29   #define PetscCallPTSCOTCH(...) \
30     do { \
31       PetscCheck(!(__VA_ARGS__), PETSC_COMM_SELF, PETSC_ERR_LIB, "Error calling PT-Scotch library"); \
32     } while (0)
33 
PTScotch_Strategy(PetscInt strategy)34 static int PTScotch_Strategy(PetscInt strategy)
35 {
36   switch (strategy) {
37   case 0:
38     return SCOTCH_STRATDEFAULT;
39   case 1:
40     return SCOTCH_STRATQUALITY;
41   case 2:
42     return SCOTCH_STRATSPEED;
43   case 3:
44     return SCOTCH_STRATBALANCE;
45   case 4:
46     return SCOTCH_STRATSAFETY;
47   case 5:
48     return SCOTCH_STRATSCALABILITY;
49   case 6:
50     return SCOTCH_STRATRECURSIVE;
51   case 7:
52     return SCOTCH_STRATREMAP;
53   default:
54     return SCOTCH_STRATDEFAULT;
55   }
56 }
57 
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[])58 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[])
59 {
60   SCOTCH_Arch  archdat;
61   SCOTCH_Graph grafdat;
62   SCOTCH_Strat stradat;
63   SCOTCH_Num   vertnbr = n;
64   SCOTCH_Num   edgenbr = xadj[n];
65   SCOTCH_Num  *velotab = vtxwgt;
66   SCOTCH_Num  *edlotab = adjwgt;
67   SCOTCH_Num   flagval = strategy;
68   double       kbalval = imbalance;
69 
70   PetscFunctionBegin;
71   if (!n) PetscFunctionReturn(PETSC_SUCCESS);
72   {
73     PetscBool flg = PETSC_TRUE;
74     PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", "-petscpartitioner_use_vertex_weights", "3.13", NULL));
75     /*
76        Cannot remove the PetscOptionsGetBool() below since the PetscOptionsDeprecatedNoObject() above is called after the non-deprecated version
77        has already been checked in PetscPartitionerSetFromOptions().
78     */
79     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_use_vertex_weight", &flg, NULL));
80     if (!flg) velotab = NULL;
81   }
82   PetscCallPTSCOTCH(SCOTCH_graphInit(&grafdat));
83   PetscCallPTSCOTCH(SCOTCH_graphBuild(&grafdat, 0, vertnbr, xadj, xadj + 1, velotab, NULL, edgenbr, adjncy, edlotab));
84   PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat));
85   PetscCallPTSCOTCH(SCOTCH_stratGraphMapBuild(&stradat, flagval, nparts, kbalval));
86   PetscCallPTSCOTCH(SCOTCH_archInit(&archdat));
87   if (tpart) {
88     PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, PetscMin(nparts, n), tpart));
89   } else {
90     PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, PetscMin(nparts, n)));
91   }
92   PetscCallPTSCOTCH(SCOTCH_graphMap(&grafdat, &archdat, &stradat, part));
93   SCOTCH_archExit(&archdat);
94   SCOTCH_stratExit(&stradat);
95   SCOTCH_graphExit(&grafdat);
96   PetscFunctionReturn(PETSC_SUCCESS);
97 }
98 
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)99 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)
100 {
101   PetscMPIInt     procglbnbr;
102   PetscMPIInt     proclocnum;
103   SCOTCH_Arch     archdat;
104   SCOTCH_Dgraph   grafdat;
105   SCOTCH_Dmapping mappdat;
106   SCOTCH_Strat    stradat;
107   SCOTCH_Num      vertlocnbr;
108   SCOTCH_Num      edgelocnbr;
109   SCOTCH_Num     *veloloctab = vtxwgt;
110   SCOTCH_Num     *edloloctab = adjwgt;
111   SCOTCH_Num      flagval    = strategy;
112   double          kbalval    = imbalance;
113 
114   PetscFunctionBegin;
115   {
116     PetscBool flg = PETSC_TRUE;
117     PetscCall(PetscOptionsDeprecatedNoObject("-petscpartititoner_ptscotch_vertex_weight", "-petscpartitioner_use_vertex_weights", "3.13", NULL));
118     /*
119        Cannot remove the PetscOptionsGetBool() below since the PetscOptionsDeprecatedNoObject() above is called after the non-deprecated version
120        has already been checked in PetscPartitionerSetFromOptions().
121     */
122     PetscCall(PetscOptionsGetBool(NULL, NULL, "-petscpartititoner_use_vertex_weight", &flg, NULL));
123     if (!flg) veloloctab = NULL;
124   }
125   PetscCallMPI(MPI_Comm_size(comm, &procglbnbr));
126   PetscCallMPI(MPI_Comm_rank(comm, &proclocnum));
127   vertlocnbr = vtxdist[proclocnum + 1] - vtxdist[proclocnum];
128   edgelocnbr = xadj[vertlocnbr];
129 
130   PetscCallPTSCOTCH(SCOTCH_dgraphInit(&grafdat, comm));
131   PetscCallPTSCOTCH(SCOTCH_dgraphBuild(&grafdat, 0, vertlocnbr, vertlocnbr, xadj, xadj + 1, veloloctab, NULL, edgelocnbr, edgelocnbr, adjncy, NULL, edloloctab));
132   PetscCallPTSCOTCH(SCOTCH_stratInit(&stradat));
133   PetscCallPTSCOTCH(SCOTCH_stratDgraphMapBuild(&stradat, flagval, procglbnbr, nparts, kbalval));
134   PetscCallPTSCOTCH(SCOTCH_archInit(&archdat));
135   if (tpart) { /* target partition weights */
136     PetscCallPTSCOTCH(SCOTCH_archCmpltw(&archdat, nparts, tpart));
137   } else {
138     PetscCallPTSCOTCH(SCOTCH_archCmplt(&archdat, nparts));
139   }
140   PetscCallPTSCOTCH(SCOTCH_dgraphMapInit(&grafdat, &mappdat, &archdat, part));
141   PetscCallPTSCOTCH(SCOTCH_dgraphMapCompute(&grafdat, &mappdat, &stradat));
142   SCOTCH_dgraphMapExit(&grafdat, &mappdat);
143   SCOTCH_archExit(&archdat);
144   SCOTCH_stratExit(&stradat);
145   SCOTCH_dgraphExit(&grafdat);
146   PetscFunctionReturn(PETSC_SUCCESS);
147 }
148 
149 #endif /* PETSC_HAVE_PTSCOTCH */
150 
151 static const char *const PTScotchStrategyList[] = {"DEFAULT", "QUALITY", "SPEED", "BALANCE", "SAFETY", "SCALABILITY", "RECURSIVE", "REMAP"};
152 
PetscPartitionerDestroy_PTScotch(PetscPartitioner part)153 static PetscErrorCode PetscPartitionerDestroy_PTScotch(PetscPartitioner part)
154 {
155   PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data;
156 
157   PetscFunctionBegin;
158   PetscCallMPI(MPI_Comm_free(&p->pcomm));
159   PetscCall(PetscFree(part->data));
160   PetscFunctionReturn(PETSC_SUCCESS);
161 }
162 
PetscPartitionerView_PTScotch_ASCII(PetscPartitioner part,PetscViewer viewer)163 static PetscErrorCode PetscPartitionerView_PTScotch_ASCII(PetscPartitioner part, PetscViewer viewer)
164 {
165   PetscPartitioner_PTScotch *p = (PetscPartitioner_PTScotch *)part->data;
166 
167   PetscFunctionBegin;
168   PetscCall(PetscViewerASCIIPushTab(viewer));
169   PetscCall(PetscViewerASCIIPrintf(viewer, "using partitioning strategy %s\n", PTScotchStrategyList[p->strategy]));
170   PetscCall(PetscViewerASCIIPrintf(viewer, "using load imbalance ratio %g\n", (double)p->imbalance));
171   PetscCall(PetscViewerASCIIPopTab(viewer));
172   PetscFunctionReturn(PETSC_SUCCESS);
173 }
174 
PetscPartitionerView_PTScotch(PetscPartitioner part,PetscViewer viewer)175 static PetscErrorCode PetscPartitionerView_PTScotch(PetscPartitioner part, PetscViewer viewer)
176 {
177   PetscBool isascii;
178 
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
181   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
182   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
183   if (isascii) PetscCall(PetscPartitionerView_PTScotch_ASCII(part, viewer));
184   PetscFunctionReturn(PETSC_SUCCESS);
185 }
186 
PetscPartitionerSetFromOptions_PTScotch(PetscPartitioner part,PetscOptionItems PetscOptionsObject)187 static PetscErrorCode PetscPartitionerSetFromOptions_PTScotch(PetscPartitioner part, PetscOptionItems PetscOptionsObject)
188 {
189   PetscPartitioner_PTScotch *p     = (PetscPartitioner_PTScotch *)part->data;
190   const char *const         *slist = PTScotchStrategyList;
191   PetscInt                   nlist = PETSC_STATIC_ARRAY_LENGTH(PTScotchStrategyList);
192   PetscBool                  flag;
193 
194   PetscFunctionBegin;
195   PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner PTScotch Options");
196   PetscCall(PetscOptionsEList("-petscpartitioner_ptscotch_strategy", "Partitioning strategy", "", slist, nlist, slist[p->strategy], &p->strategy, &flag));
197   PetscCall(PetscOptionsReal("-petscpartitioner_ptscotch_imbalance", "Load imbalance ratio", "", p->imbalance, &p->imbalance, &flag));
198   PetscOptionsHeadEnd();
199   PetscFunctionReturn(PETSC_SUCCESS);
200 }
201 
PetscPartitionerPartition_PTScotch(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * partition)202 static PetscErrorCode PetscPartitionerPartition_PTScotch(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
203 {
204 #if defined(PETSC_HAVE_PTSCOTCH)
205   MPI_Comm    comm;
206   PetscInt    nvtxs = numVertices; /* The number of vertices in full graph */
207   PetscInt   *vtxdist;             /* Distribution of vertices across processes */
208   PetscInt   *xadj   = start;      /* Start of edge list for each vertex */
209   PetscInt   *adjncy = adjacency;  /* Edge lists for all vertices */
210   PetscInt   *vwgt   = NULL;       /* Vertex weights */
211   PetscInt   *adjwgt = NULL;       /* Edge weights */
212   PetscInt    v, i, *assignment, *points;
213   PetscMPIInt size, rank, p;
214   PetscBool   hasempty = PETSC_FALSE;
215   PetscInt   *tpwgts   = NULL;
216 
217   PetscFunctionBegin;
218   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
219   PetscCallMPI(MPI_Comm_size(comm, &size));
220   PetscCallMPI(MPI_Comm_rank(comm, &rank));
221   PetscCall(PetscMalloc2(size + 1, &vtxdist, PetscMax(nvtxs, 1), &assignment));
222   /* Calculate vertex distribution */
223   vtxdist[0] = 0;
224   PetscCallMPI(MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm));
225   for (p = 2; p <= size; ++p) {
226     hasempty = (PetscBool)(hasempty || !vtxdist[p - 1] || !vtxdist[p]);
227     vtxdist[p] += vtxdist[p - 1];
228   }
229   /* null graph */
230   if (vtxdist[size] == 0) {
231     PetscCall(PetscFree2(vtxdist, assignment));
232     PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition));
233     PetscFunctionReturn(PETSC_SUCCESS);
234   }
235 
236   /* Calculate vertex weights */
237   if (vertSection) {
238     PetscCall(PetscMalloc1(nvtxs, &vwgt));
239     for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v]));
240   }
241   // Weight edges
242   if (edgeSection) {
243     PetscCall(PetscMalloc1(xadj[nvtxs], &adjwgt));
244     for (PetscInt e = 0; e < xadj[nvtxs]; ++e) PetscCall(PetscSectionGetDof(edgeSection, e, &adjwgt[e]));
245   }
246 
247   /* Calculate partition weights */
248   if (targetSection) {
249     PetscInt sumw;
250 
251     PetscCall(PetscCalloc1(nparts, &tpwgts));
252     for (p = 0, sumw = 0; p < nparts; ++p) {
253       PetscCall(PetscSectionGetDof(targetSection, p, &tpwgts[p]));
254       sumw += tpwgts[p];
255     }
256     if (!sumw) PetscCall(PetscFree(tpwgts));
257   }
258 
259   {
260     PetscPartitioner_PTScotch *pts   = (PetscPartitioner_PTScotch *)part->data;
261     int                        strat = PTScotch_Strategy(pts->strategy);
262     double                     imbal = (double)pts->imbalance;
263 
264     for (p = 0; !vtxdist[p + 1] && p < size; ++p);
265     if (vtxdist[p + 1] == vtxdist[size]) {
266       if (rank == p) PetscCall(PTScotch_PartGraph_Seq(strat, imbal, nvtxs, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment));
267     } else {
268       MPI_Comm pcomm = pts->pcomm;
269 
270       if (hasempty) {
271         PetscInt cnt;
272 
273         PetscCallMPI(MPI_Comm_split(pts->pcomm, !!nvtxs, rank, &pcomm));
274         for (p = 0, cnt = 0; p < size; p++) {
275           if (vtxdist[p + 1] != vtxdist[p]) {
276             vtxdist[cnt + 1] = vtxdist[p + 1];
277             cnt++;
278           }
279         }
280       }
281       if (nvtxs) PetscCall(PTScotch_PartGraph_MPI(strat, imbal, vtxdist, xadj, adjncy, vwgt, adjwgt, nparts, tpwgts, assignment, pcomm));
282       if (hasempty) PetscCallMPI(MPI_Comm_free(&pcomm));
283     }
284   }
285   PetscCall(PetscFree(vwgt));
286   PetscCall(PetscFree(adjwgt));
287   PetscCall(PetscFree(tpwgts));
288 
289   /* Convert to PetscSection+IS */
290   for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1));
291   PetscCall(PetscMalloc1(nvtxs, &points));
292   for (p = 0, i = 0; p < nparts; ++p) {
293     for (v = 0; v < nvtxs; ++v) {
294       if (assignment[v] == p) points[i++] = v;
295     }
296   }
297   PetscCheck(i == nvtxs, comm, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs);
298   PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition));
299 
300   PetscCall(PetscFree2(vtxdist, assignment));
301   PetscFunctionReturn(PETSC_SUCCESS);
302 #else
303   SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-ptscotch.");
304 #endif
305 }
306 
PetscPartitionerInitialize_PTScotch(PetscPartitioner part)307 static PetscErrorCode PetscPartitionerInitialize_PTScotch(PetscPartitioner part)
308 {
309   PetscFunctionBegin;
310   part->noGraph             = PETSC_FALSE;
311   part->ops->view           = PetscPartitionerView_PTScotch;
312   part->ops->destroy        = PetscPartitionerDestroy_PTScotch;
313   part->ops->partition      = PetscPartitionerPartition_PTScotch;
314   part->ops->setfromoptions = PetscPartitionerSetFromOptions_PTScotch;
315   PetscFunctionReturn(PETSC_SUCCESS);
316 }
317 
318 /*MC
319   PETSCPARTITIONERPTSCOTCH = "ptscotch" - A PetscPartitioner object using the PT-Scotch library
320 
321   Level: intermediate
322 
323   Options Database Keys:
324 +  -petscpartitioner_ptscotch_strategy <string> - PT-Scotch strategy. Choose one of default quality speed balance safety scalability recursive remap
325 -  -petscpartitioner_ptscotch_imbalance <val> - Load imbalance ratio
326 
327   Notes: when the graph is on a single process, this partitioner actually uses Scotch and not PT-Scotch
328 
329 .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
330 M*/
331 
PetscPartitionerCreate_PTScotch(PetscPartitioner part)332 PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_PTScotch(PetscPartitioner part)
333 {
334   PetscPartitioner_PTScotch *p;
335 
336   PetscFunctionBegin;
337   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
338   PetscCall(PetscNew(&p));
339   part->data = p;
340 
341   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)part), &p->pcomm));
342   p->strategy  = 0;
343   p->imbalance = 0.01;
344 
345   PetscCall(PetscPartitionerInitialize_PTScotch(part));
346   PetscCall(PetscCitationsRegister(PTScotchPartitionerCitation, &PTScotchPartitionerCite));
347   PetscFunctionReturn(PETSC_SUCCESS);
348 }
349