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