xref: /petsc/src/mat/graphops/partition/impls/party/party.c (revision 8838bf166a5dcd11a3e89db66de7a74944683c57)
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 
MatPartitioningApply_Party(MatPartitioning part,IS * partitioning)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 
MatPartitioningView_Party(MatPartitioning part,PetscViewer viewer)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 /*@
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   Developer Note:
173   Should be `MatPartitioningPartySetGlobalType()` and all uses of method should be changed to type
174 
175 .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetLocal()`
176 @*/
MatPartitioningPartySetGlobal(MatPartitioning part,const char * global)177 PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part, const char *global)
178 {
179   PetscFunctionBegin;
180   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
181   PetscTryMethod(part, "MatPartitioningPartySetGlobal_C", (MatPartitioning, const char *), (part, global));
182   PetscFunctionReturn(PETSC_SUCCESS);
183 }
184 
MatPartitioningPartySetGlobal_Party(MatPartitioning part,const char * global)185 static PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part, const char *global)
186 {
187   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
188 
189   PetscFunctionBegin;
190   PetscCall(PetscStrncpy(party->global, global, 15));
191   PetscFunctionReturn(PETSC_SUCCESS);
192 }
193 
194 /*@
195   MatPartitioningPartySetLocal - Set local method used by the Party partitioner.
196 
197   Collective
198 
199   Input Parameters:
200 + part  - the partitioning context
201 - local - a string representing the method
202 
203   Options Database Key:
204 . -mat_partitioning_party_local <method> - the local method
205 
206   Level: advanced
207 
208   Note:
209   The method may be one of `MP_PARTY_HELPFUL_SETS`, `MP_PARTY_KERNIGHAN_LIN`, or
210   `MP_PARTY_NONE`. Check the Party Library Users Manual for details.
211 
212   Developer Note:
213   Should be `MatPartitioningPartySetLocalType()` and all uses of method should be changed to type
214 
215 .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetGlobal()`
216 @*/
MatPartitioningPartySetLocal(MatPartitioning part,const char * local)217 PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part, const char *local)
218 {
219   PetscFunctionBegin;
220   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
221   PetscTryMethod(part, "MatPartitioningPartySetLocal_C", (MatPartitioning, const char *), (part, local));
222   PetscFunctionReturn(PETSC_SUCCESS);
223 }
224 
MatPartitioningPartySetLocal_Party(MatPartitioning part,const char * local)225 static PetscErrorCode MatPartitioningPartySetLocal_Party(MatPartitioning part, const char *local)
226 {
227   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
228 
229   PetscFunctionBegin;
230   PetscCall(PetscStrncpy(party->local, local, 15));
231   PetscFunctionReturn(PETSC_SUCCESS);
232 }
233 
234 /*@
235   MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
236   Party partitioner.
237 
238   Collective
239 
240   Input Parameters:
241 + part  - the partitioning context
242 - level - the coarse level in range [0.0,1.0]
243 
244   Options Database Key:
245 . -mat_partitioning_party_coarse <l> - Coarse level
246 
247   Level: advanced
248 
249 .seealso: `MATPARTITIONINGPARTY`
250 @*/
MatPartitioningPartySetCoarseLevel(MatPartitioning part,PetscReal level)251 PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part, PetscReal level)
252 {
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
255   PetscValidLogicalCollectiveReal(part, level, 2);
256   PetscTryMethod(part, "MatPartitioningPartySetCoarseLevel_C", (MatPartitioning, PetscReal), (part, level));
257   PetscFunctionReturn(PETSC_SUCCESS);
258 }
259 
MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part,PetscReal level)260 static PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part, PetscReal level)
261 {
262   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
263 
264   PetscFunctionBegin;
265   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]");
266   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
267   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
268   PetscFunctionReturn(PETSC_SUCCESS);
269 }
270 
271 /*@
272   MatPartitioningPartySetMatchOptimization - Activate matching optimization for
273   graph reduction.
274 
275   Collective
276 
277   Input Parameters:
278 + part - the partitioning context
279 - opt  - boolean flag
280 
281   Options Database Key:
282 . -mat_partitioning_party_match_optimization - Matching optimization on/off
283 
284   Level: advanced
285 
286 .seealso: `MATPARTITIONINGPARTY`
287 @*/
MatPartitioningPartySetMatchOptimization(MatPartitioning part,PetscBool opt)288 PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part, PetscBool opt)
289 {
290   PetscFunctionBegin;
291   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
292   PetscValidLogicalCollectiveBool(part, opt, 2);
293   PetscTryMethod(part, "MatPartitioningPartySetMatchOptimization_C", (MatPartitioning, PetscBool), (part, opt));
294   PetscFunctionReturn(PETSC_SUCCESS);
295 }
296 
MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part,PetscBool opt)297 static PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part, PetscBool opt)
298 {
299   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
300 
301   PetscFunctionBegin;
302   party->redo = opt;
303   PetscFunctionReturn(PETSC_SUCCESS);
304 }
305 
306 /*@
307   MatPartitioningPartySetBipart - Activate or deactivate recursive bisection in the Party partitioner
308 
309   Collective
310 
311   Input Parameters:
312 + part - the partitioning context
313 - bp   - boolean flag
314 
315   Options Database Key:
316 . -mat_partitioning_party_bipart - Bipartitioning option on/off
317 
318   Level: advanced
319 
320 .seealso: `MATPARTITIONINGPARTY`
321 @*/
MatPartitioningPartySetBipart(MatPartitioning part,PetscBool bp)322 PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part, PetscBool bp)
323 {
324   PetscFunctionBegin;
325   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
326   PetscValidLogicalCollectiveBool(part, bp, 2);
327   PetscTryMethod(part, "MatPartitioningPartySetBipart_C", (MatPartitioning, PetscBool), (part, bp));
328   PetscFunctionReturn(PETSC_SUCCESS);
329 }
330 
MatPartitioningPartySetBipart_Party(MatPartitioning part,PetscBool bp)331 static PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part, PetscBool bp)
332 {
333   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
334 
335   PetscFunctionBegin;
336   party->recursive = bp;
337   PetscFunctionReturn(PETSC_SUCCESS);
338 }
339 
MatPartitioningSetFromOptions_Party(MatPartitioning part,PetscOptionItems PetscOptionsObject)340 static PetscErrorCode MatPartitioningSetFromOptions_Party(MatPartitioning part, PetscOptionItems PetscOptionsObject)
341 {
342   PetscBool              flag;
343   char                   value[256];
344   PetscReal              r;
345   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
346 
347   PetscFunctionBegin;
348   PetscOptionsHeadBegin(PetscOptionsObject, "Set Party partitioning options");
349   PetscCall(PetscOptionsString("-mat_partitioning_party_global", "Global method", "MatPartitioningPartySetGlobal", party->global, value, sizeof(value), &flag));
350   if (flag) PetscCall(MatPartitioningPartySetGlobal(part, value));
351   PetscCall(PetscOptionsString("-mat_partitioning_party_local", "Local method", "MatPartitioningPartySetLocal", party->local, value, sizeof(value), &flag));
352   if (flag) PetscCall(MatPartitioningPartySetLocal(part, value));
353   PetscCall(PetscOptionsReal("-mat_partitioning_party_coarse", "Coarse level", "MatPartitioningPartySetCoarseLevel", 0.0, &r, &flag));
354   if (flag) PetscCall(MatPartitioningPartySetCoarseLevel(part, r));
355   PetscCall(PetscOptionsBool("-mat_partitioning_party_match_optimization", "Matching optimization on/off", "MatPartitioningPartySetMatchOptimization", party->redo, &party->redo, NULL));
356   PetscCall(PetscOptionsBool("-mat_partitioning_party_bipart", "Bipartitioning on/off", "MatPartitioningPartySetBipart", party->recursive, &party->recursive, NULL));
357   PetscCall(PetscOptionsBool("-mat_partitioning_party_verbose", "Show library output", "", party->verbose, &party->verbose, NULL));
358   PetscOptionsHeadEnd();
359   PetscFunctionReturn(PETSC_SUCCESS);
360 }
361 
MatPartitioningDestroy_Party(MatPartitioning part)362 static PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
363 {
364   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
365 
366   PetscFunctionBegin;
367   PetscCall(PetscFree(party));
368   /* clear composed functions */
369   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", NULL));
370   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", NULL));
371   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", NULL));
372   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", NULL));
373   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", NULL));
374   PetscFunctionReturn(PETSC_SUCCESS);
375 }
376 
377 /*MC
378    MATPARTITIONINGPARTY - Creates a partitioning context via the external package Party <
379    http://wwwcs.upb.de/fachbereich/AG/monien/RESEARCH/PART/party.htm>.
380 
381    Level: beginner
382 
383    Note:
384    Does not support the `MatPartitioningSetUseEdgeWeights()` option
385 
386 .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningPartySetGlobal()`, `MatPartitioningPartySetLocal()`,
387           `MatPartitioningPartySetCoarseLevel()`, `MatPartitioningPartySetMatchOptimization()`, `MatPartitioningPartySetBipart()`
388 M*/
389 
MatPartitioningCreate_Party(MatPartitioning part)390 PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning part)
391 {
392   MatPartitioning_Party *party;
393 
394   PetscFunctionBegin;
395   PetscCall(PetscNew(&party));
396   part->data = (void *)party;
397 
398   PetscCall(PetscStrncpy(party->global, "gcf,gbf", sizeof(party->global)));
399   PetscCall(PetscStrncpy(party->local, "kl", sizeof(party->local)));
400 
401   party->redm         = PETSC_TRUE;
402   party->redo         = PETSC_TRUE;
403   party->recursive    = PETSC_TRUE;
404   party->verbose      = PETSC_FALSE;
405   party->nbvtxcoarsed = 200;
406 
407   part->ops->apply          = MatPartitioningApply_Party;
408   part->ops->view           = MatPartitioningView_Party;
409   part->ops->destroy        = MatPartitioningDestroy_Party;
410   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;
411 
412   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", MatPartitioningPartySetGlobal_Party));
413   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", MatPartitioningPartySetLocal_Party));
414   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", MatPartitioningPartySetCoarseLevel_Party));
415   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", MatPartitioningPartySetMatchOptimization_Party));
416   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", MatPartitioningPartySetBipart_Party));
417   PetscFunctionReturn(PETSC_SUCCESS);
418 }
419