xref: /petsc/src/mat/graphops/partition/impls/party/party.c (revision d9acb416d05abeed0a33bde3a81aeb2ea0364f6a)
1 #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
2 
3 #if defined(PETSC_HAVE_UNISTD_H)
4   #include <unistd.h>
5 #endif
6 
7 EXTERN_C_BEGIN
8 #include <party_lib.h>
9 EXTERN_C_END
10 
11 typedef struct {
12   PetscBool redm;
13   PetscBool redo;
14   PetscBool recursive;
15   PetscBool verbose;
16   char      global[15];   /* global method */
17   char      local[15];    /* local method */
18   PetscInt  nbvtxcoarsed; /* number of vertices for the coarse graph */
19 } MatPartitioning_Party;
20 
21 #define SIZE_LOG 10000 /* size of buffer for mesg_log */
22 
23 static PetscErrorCode MatPartitioningApply_Party(MatPartitioning part, IS *partitioning)
24 {
25   int                    perr;
26   PetscInt               i, *parttab, *locals, nb_locals, M, N;
27   PetscMPIInt            size, rank;
28   Mat                    mat = part->adj, matAdj, matSeq, *A;
29   Mat_MPIAdj            *adj;
30   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
31   PetscBool              flg;
32   IS                     isrow, iscol;
33   int                    n, *edge_p, *edge, *vertex_w, p, *part_party, cutsize, redl, rec;
34   const char            *redm, *redo;
35   char                  *mesg_log;
36 #if defined(PETSC_HAVE_UNISTD_H)
37   int fd_stdout, fd_pipe[2], count;
38 #endif
39 
40   PetscFunctionBegin;
41   PetscCheck(!part->use_edge_weights, PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Party does not support edge weights");
42   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
43   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
44   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
45   if (size > 1) {
46     if (flg) {
47       PetscCall(MatMPIAdjToSeq(mat, &matSeq));
48     } else {
49       PetscCall(PetscInfo(part, "Converting distributed matrix to sequential: this could be a performance loss\n"));
50       PetscCall(MatGetSize(mat, &M, &N));
51       PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &isrow));
52       PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol));
53       PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &A));
54       PetscCall(ISDestroy(&isrow));
55       PetscCall(ISDestroy(&iscol));
56       matSeq = *A;
57       PetscCall(PetscFree(A));
58     }
59   } else {
60     PetscCall(PetscObjectReference((PetscObject)mat));
61     matSeq = mat;
62   }
63 
64   if (!flg) { /* convert regular matrix to MPIADJ */
65     PetscCall(MatConvert(matSeq, MATMPIADJ, MAT_INITIAL_MATRIX, &matAdj));
66   } else {
67     PetscCall(PetscObjectReference((PetscObject)matSeq));
68     matAdj = matSeq;
69   }
70 
71   adj = (Mat_MPIAdj *)matAdj->data; /* finally adj contains adjacency graph */
72 
73   /* arguments for Party library */
74   n        = mat->rmap->N;             /* number of vertices in full graph */
75   edge_p   = adj->i;                   /* start of edge list for each vertex */
76   edge     = adj->j;                   /* edge list data */
77   vertex_w = part->vertex_weights;     /* weights for all vertices */
78   p        = part->n;                  /* number of parts to create */
79   redl     = party->nbvtxcoarsed;      /* how many vertices to coarsen down to? */
80   rec      = party->recursive ? 1 : 0; /* recursive bisection */
81   redm     = party->redm ? "lam" : ""; /* matching method */
82   redo     = party->redo ? "w3" : "";  /* matching optimization method */
83 
84   PetscCall(PetscMalloc1(mat->rmap->N, &part_party));
85 
86   /* redirect output to buffer */
87 #if defined(PETSC_HAVE_UNISTD_H)
88   fd_stdout = dup(1);
89   PetscCheck(!pipe(fd_pipe), PETSC_COMM_SELF, PETSC_ERR_SYS, "Could not open pipe");
90   close(1);
91   dup2(fd_pipe[1], 1);
92   PetscCall(PetscMalloc1(SIZE_LOG, &mesg_log));
93 #endif
94 
95   /* library call */
96   party_lib_times_start();
97   perr = party_lib(n, vertex_w, NULL, NULL, NULL, edge_p, edge, NULL, p, part_party, &cutsize, redl, (char *)redm, (char *)redo, party->global, party->local, rec, 1);
98 
99   party_lib_times_output(1);
100   part_info(n, vertex_w, edge_p, edge, NULL, p, part_party, 1);
101 
102 #if defined(PETSC_HAVE_UNISTD_H)
103   PetscCall(PetscFFlush(stdout));
104   count = read(fd_pipe[0], mesg_log, (SIZE_LOG - 1) * sizeof(char));
105   if (count < 0) count = 0;
106   mesg_log[count] = 0;
107   close(1);
108   dup2(fd_stdout, 1);
109   close(fd_stdout);
110   close(fd_pipe[0]);
111   close(fd_pipe[1]);
112   if (party->verbose) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "%s", mesg_log));
113   PetscCall(PetscFree(mesg_log));
114 #endif
115   PetscCheck(!perr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Party failed");
116 
117   PetscCall(PetscMalloc1(mat->rmap->N, &parttab));
118   for (i = 0; i < mat->rmap->N; i++) parttab[i] = part_party[i];
119 
120   /* creation of the index set */
121   nb_locals = mat->rmap->n;
122   locals    = parttab + mat->rmap->rstart;
123 
124   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), nb_locals, locals, PETSC_COPY_VALUES, partitioning));
125 
126   /* clean up */
127   PetscCall(PetscFree(parttab));
128   PetscCall(PetscFree(part_party));
129   PetscCall(MatDestroy(&matSeq));
130   PetscCall(MatDestroy(&matAdj));
131   PetscFunctionReturn(PETSC_SUCCESS);
132 }
133 
134 static PetscErrorCode MatPartitioningView_Party(MatPartitioning part, PetscViewer viewer)
135 {
136   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
137   PetscBool              isascii;
138 
139   PetscFunctionBegin;
140   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
141   if (isascii) {
142     PetscCall(PetscViewerASCIIPrintf(viewer, "  Global method: %s\n", party->global));
143     PetscCall(PetscViewerASCIIPrintf(viewer, "  Local method: %s\n", party->local));
144     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of vertices for the coarse graph: %d\n", party->nbvtxcoarsed));
145     if (party->redm) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using matching method for graph reduction\n"));
146     if (party->redo) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using matching optimization\n"));
147     if (party->recursive) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using recursive bipartitioning\n"));
148   }
149   PetscFunctionReturn(PETSC_SUCCESS);
150 }
151 
152 /*@C
153   MatPartitioningPartySetGlobal - Set global method for Party partitioner.
154 
155   Collective
156 
157   Input Parameters:
158 + part   - the partitioning context
159 - global - a string representing the method
160 
161   Options Database Key:
162 . -mat_partitioning_party_global <method> - the global method
163 
164   Level: advanced
165 
166   Note:
167   The method may be one of `MP_PARTY_OPT`, `MP_PARTY_LIN`, `MP_PARTY_SCA`,
168   `MP_PARTY_RAN`, `MP_PARTY_GBF`, `MP_PARTY_GCF`, `MP_PARTY_BUB` or `MP_PARTY_DEF`, or
169   alternatively a string describing the method. Two or more methods can be
170   combined like "gbf,gcf". Check the Party Library Users Manual for details.
171 
172 .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetLocal()`
173 @*/
174 PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part, const char *global)
175 {
176   PetscFunctionBegin;
177   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
178   PetscTryMethod(part, "MatPartitioningPartySetGlobal_C", (MatPartitioning, const char *), (part, global));
179   PetscFunctionReturn(PETSC_SUCCESS);
180 }
181 
182 static PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part, const char *global)
183 {
184   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
185 
186   PetscFunctionBegin;
187   PetscCall(PetscStrncpy(party->global, global, 15));
188   PetscFunctionReturn(PETSC_SUCCESS);
189 }
190 
191 /*@C
192   MatPartitioningPartySetLocal - Set local method used by the Party partitioner.
193 
194   Collective
195 
196   Input Parameters:
197 + part  - the partitioning context
198 - local - a string representing the method
199 
200   Options Database Key:
201 . -mat_partitioning_party_local <method> - the local method
202 
203   Level: advanced
204 
205   Note:
206   The method may be one of `MP_PARTY_HELPFUL_SETS`, `MP_PARTY_KERNIGHAN_LIN`, or
207   `MP_PARTY_NONE`. Check the Party Library Users Manual for details.
208 
209 .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetGlobal()`
210 @*/
211 PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part, const char *local)
212 {
213   PetscFunctionBegin;
214   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
215   PetscTryMethod(part, "MatPartitioningPartySetLocal_C", (MatPartitioning, const char *), (part, local));
216   PetscFunctionReturn(PETSC_SUCCESS);
217 }
218 
219 static PetscErrorCode MatPartitioningPartySetLocal_Party(MatPartitioning part, const char *local)
220 
221 {
222   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
223 
224   PetscFunctionBegin;
225   PetscCall(PetscStrncpy(party->local, local, 15));
226   PetscFunctionReturn(PETSC_SUCCESS);
227 }
228 
229 /*@
230   MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
231   Party partitioner.
232 
233   Collective
234 
235   Input Parameters:
236 + part  - the partitioning context
237 - level - the coarse level in range [0.0,1.0]
238 
239   Options Database Key:
240 . -mat_partitioning_party_coarse <l> - Coarse level
241 
242   Level: advanced
243 
244 .seealso: `MATPARTITIONINGPARTY`
245 @*/
246 PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part, PetscReal level)
247 {
248   PetscFunctionBegin;
249   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
250   PetscValidLogicalCollectiveReal(part, level, 2);
251   PetscTryMethod(part, "MatPartitioningPartySetCoarseLevel_C", (MatPartitioning, PetscReal), (part, level));
252   PetscFunctionReturn(PETSC_SUCCESS);
253 }
254 
255 static PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part, PetscReal level)
256 {
257   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
258 
259   PetscFunctionBegin;
260   PetscCheck(level >= 0.0 && level <= 1.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Party: level of coarsening out of range [0.0-1.0]");
261   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
262   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
263   PetscFunctionReturn(PETSC_SUCCESS);
264 }
265 
266 /*@
267   MatPartitioningPartySetMatchOptimization - Activate matching optimization for
268   graph reduction.
269 
270   Collective
271 
272   Input Parameters:
273 + part - the partitioning context
274 - opt  - boolean flag
275 
276   Options Database Key:
277 . -mat_partitioning_party_match_optimization - Matching optimization on/off
278 
279   Level: advanced
280 
281 .seealso: `MATPARTITIONINGPARTY`
282 @*/
283 PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part, PetscBool opt)
284 {
285   PetscFunctionBegin;
286   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
287   PetscValidLogicalCollectiveBool(part, opt, 2);
288   PetscTryMethod(part, "MatPartitioningPartySetMatchOptimization_C", (MatPartitioning, PetscBool), (part, opt));
289   PetscFunctionReturn(PETSC_SUCCESS);
290 }
291 
292 static PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part, PetscBool opt)
293 {
294   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
295 
296   PetscFunctionBegin;
297   party->redo = opt;
298   PetscFunctionReturn(PETSC_SUCCESS);
299 }
300 
301 /*@
302   MatPartitioningPartySetBipart - Activate or deactivate recursive bisection in the Party partitioner
303 
304   Collective
305 
306   Input Parameters:
307 + part - the partitioning context
308 - bp   - boolean flag
309 
310   Options Database Key:
311 . -mat_partitioning_party_bipart - Bipartitioning option on/off
312 
313   Level: advanced
314 
315 .seealso: `MATPARTITIONINGPARTY`
316 @*/
317 PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part, PetscBool bp)
318 {
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
321   PetscValidLogicalCollectiveBool(part, bp, 2);
322   PetscTryMethod(part, "MatPartitioningPartySetBipart_C", (MatPartitioning, PetscBool), (part, bp));
323   PetscFunctionReturn(PETSC_SUCCESS);
324 }
325 
326 static PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part, PetscBool bp)
327 {
328   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
329 
330   PetscFunctionBegin;
331   party->recursive = bp;
332   PetscFunctionReturn(PETSC_SUCCESS);
333 }
334 
335 static PetscErrorCode MatPartitioningSetFromOptions_Party(MatPartitioning part, PetscOptionItems *PetscOptionsObject)
336 {
337   PetscBool              flag;
338   char                   value[256];
339   PetscReal              r;
340   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
341 
342   PetscFunctionBegin;
343   PetscOptionsHeadBegin(PetscOptionsObject, "Set Party partitioning options");
344   PetscCall(PetscOptionsString("-mat_partitioning_party_global", "Global method", "MatPartitioningPartySetGlobal", party->global, value, sizeof(value), &flag));
345   if (flag) PetscCall(MatPartitioningPartySetGlobal(part, value));
346   PetscCall(PetscOptionsString("-mat_partitioning_party_local", "Local method", "MatPartitioningPartySetLocal", party->local, value, sizeof(value), &flag));
347   if (flag) PetscCall(MatPartitioningPartySetLocal(part, value));
348   PetscCall(PetscOptionsReal("-mat_partitioning_party_coarse", "Coarse level", "MatPartitioningPartySetCoarseLevel", 0.0, &r, &flag));
349   if (flag) PetscCall(MatPartitioningPartySetCoarseLevel(part, r));
350   PetscCall(PetscOptionsBool("-mat_partitioning_party_match_optimization", "Matching optimization on/off", "MatPartitioningPartySetMatchOptimization", party->redo, &party->redo, NULL));
351   PetscCall(PetscOptionsBool("-mat_partitioning_party_bipart", "Bipartitioning on/off", "MatPartitioningPartySetBipart", party->recursive, &party->recursive, NULL));
352   PetscCall(PetscOptionsBool("-mat_partitioning_party_verbose", "Show library output", "", party->verbose, &party->verbose, NULL));
353   PetscOptionsHeadEnd();
354   PetscFunctionReturn(PETSC_SUCCESS);
355 }
356 
357 static PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
358 {
359   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
360 
361   PetscFunctionBegin;
362   PetscCall(PetscFree(party));
363   /* clear composed functions */
364   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", NULL));
365   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", NULL));
366   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", NULL));
367   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", NULL));
368   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", NULL));
369   PetscFunctionReturn(PETSC_SUCCESS);
370 }
371 
372 /*MC
373    MATPARTITIONINGPARTY - Creates a partitioning context via the external package Party <
374    http://wwwcs.upb.de/fachbereich/AG/monien/RESEARCH/PART/party.htm>.
375 
376    Level: beginner
377 
378    Note:
379    Does not support the `MatPartitioningSetUseEdgeWeights()` option
380 
381 .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningPartySetGlobal()`, `MatPartitioningPartySetLocal()`,
382           `MatPartitioningPartySetCoarseLevel()`, `MatPartitioningPartySetMatchOptimization()`, `MatPartitioningPartySetBipart()`
383 M*/
384 
385 PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning part)
386 {
387   MatPartitioning_Party *party;
388 
389   PetscFunctionBegin;
390   PetscCall(PetscNew(&party));
391   part->data = (void *)party;
392 
393   PetscCall(PetscStrncpy(party->global, "gcf,gbf", sizeof(party->global)));
394   PetscCall(PetscStrncpy(party->local, "kl", sizeof(party->local)));
395 
396   party->redm         = PETSC_TRUE;
397   party->redo         = PETSC_TRUE;
398   party->recursive    = PETSC_TRUE;
399   party->verbose      = PETSC_FALSE;
400   party->nbvtxcoarsed = 200;
401 
402   part->ops->apply          = MatPartitioningApply_Party;
403   part->ops->view           = MatPartitioningView_Party;
404   part->ops->destroy        = MatPartitioningDestroy_Party;
405   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;
406 
407   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", MatPartitioningPartySetGlobal_Party));
408   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", MatPartitioningPartySetLocal_Party));
409   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", MatPartitioningPartySetCoarseLevel_Party));
410   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", MatPartitioningPartySetMatchOptimization_Party));
411   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", MatPartitioningPartySetBipart_Party));
412   PetscFunctionReturn(PETSC_SUCCESS);
413 }
414