xref: /petsc/src/mat/graphops/partition/impls/party/party.c (revision 613ce9fe8f5e2bcdf7c72d914e9769b5d960fb4c)
18be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
28be712e4SBarry Smith 
38be712e4SBarry Smith #if defined(PETSC_HAVE_UNISTD_H)
48be712e4SBarry Smith   #include <unistd.h>
58be712e4SBarry Smith #endif
68be712e4SBarry Smith 
78be712e4SBarry Smith EXTERN_C_BEGIN
88be712e4SBarry Smith #include <party_lib.h>
98be712e4SBarry Smith EXTERN_C_END
108be712e4SBarry Smith 
118be712e4SBarry Smith typedef struct {
128be712e4SBarry Smith   PetscBool redm;
138be712e4SBarry Smith   PetscBool redo;
148be712e4SBarry Smith   PetscBool recursive;
158be712e4SBarry Smith   PetscBool verbose;
168be712e4SBarry Smith   char      global[15];   /* global method */
178be712e4SBarry Smith   char      local[15];    /* local method */
188be712e4SBarry Smith   PetscInt  nbvtxcoarsed; /* number of vertices for the coarse graph */
198be712e4SBarry Smith } MatPartitioning_Party;
208be712e4SBarry Smith 
218be712e4SBarry Smith #define SIZE_LOG 10000 /* size of buffer for mesg_log */
228be712e4SBarry Smith 
238be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Party(MatPartitioning part, IS *partitioning)
248be712e4SBarry Smith {
258be712e4SBarry Smith   int                    perr;
268be712e4SBarry Smith   PetscInt               i, *parttab, *locals, nb_locals, M, N;
278be712e4SBarry Smith   PetscMPIInt            size, rank;
288be712e4SBarry Smith   Mat                    mat = part->adj, matAdj, matSeq, *A;
298be712e4SBarry Smith   Mat_MPIAdj            *adj;
308be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
318be712e4SBarry Smith   PetscBool              flg;
328be712e4SBarry Smith   IS                     isrow, iscol;
338be712e4SBarry Smith   int                    n, *edge_p, *edge, *vertex_w, p, *part_party, cutsize, redl, rec;
348be712e4SBarry Smith   const char            *redm, *redo;
358be712e4SBarry Smith   char                  *mesg_log;
368be712e4SBarry Smith #if defined(PETSC_HAVE_UNISTD_H)
378be712e4SBarry Smith   int fd_stdout, fd_pipe[2], count;
388be712e4SBarry Smith #endif
398be712e4SBarry Smith 
408be712e4SBarry Smith   PetscFunctionBegin;
418be712e4SBarry Smith   PetscCheck(!part->use_edge_weights, PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Party does not support edge weights");
428be712e4SBarry Smith   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)mat), &size));
438be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)mat), &rank));
448be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
458be712e4SBarry Smith   if (size > 1) {
468be712e4SBarry Smith     if (flg) {
478be712e4SBarry Smith       PetscCall(MatMPIAdjToSeq(mat, &matSeq));
488be712e4SBarry Smith     } else {
498be712e4SBarry Smith       PetscCall(PetscInfo(part, "Converting distributed matrix to sequential: this could be a performance loss\n"));
508be712e4SBarry Smith       PetscCall(MatGetSize(mat, &M, &N));
518be712e4SBarry Smith       PetscCall(ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &isrow));
528be712e4SBarry Smith       PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol));
538be712e4SBarry Smith       PetscCall(MatCreateSubMatrices(mat, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &A));
548be712e4SBarry Smith       PetscCall(ISDestroy(&isrow));
558be712e4SBarry Smith       PetscCall(ISDestroy(&iscol));
568be712e4SBarry Smith       matSeq = *A;
578be712e4SBarry Smith       PetscCall(PetscFree(A));
588be712e4SBarry Smith     }
598be712e4SBarry Smith   } else {
608be712e4SBarry Smith     PetscCall(PetscObjectReference((PetscObject)mat));
618be712e4SBarry Smith     matSeq = mat;
628be712e4SBarry Smith   }
638be712e4SBarry Smith 
648be712e4SBarry Smith   if (!flg) { /* convert regular matrix to MPIADJ */
658be712e4SBarry Smith     PetscCall(MatConvert(matSeq, MATMPIADJ, MAT_INITIAL_MATRIX, &matAdj));
668be712e4SBarry Smith   } else {
678be712e4SBarry Smith     PetscCall(PetscObjectReference((PetscObject)matSeq));
688be712e4SBarry Smith     matAdj = matSeq;
698be712e4SBarry Smith   }
708be712e4SBarry Smith 
718be712e4SBarry Smith   adj = (Mat_MPIAdj *)matAdj->data; /* finally adj contains adjacency graph */
728be712e4SBarry Smith 
738be712e4SBarry Smith   /* arguments for Party library */
748be712e4SBarry Smith   n        = mat->rmap->N;             /* number of vertices in full graph */
758be712e4SBarry Smith   edge_p   = adj->i;                   /* start of edge list for each vertex */
768be712e4SBarry Smith   edge     = adj->j;                   /* edge list data */
778be712e4SBarry Smith   vertex_w = part->vertex_weights;     /* weights for all vertices */
788be712e4SBarry Smith   p        = part->n;                  /* number of parts to create */
798be712e4SBarry Smith   redl     = party->nbvtxcoarsed;      /* how many vertices to coarsen down to? */
808be712e4SBarry Smith   rec      = party->recursive ? 1 : 0; /* recursive bisection */
818be712e4SBarry Smith   redm     = party->redm ? "lam" : ""; /* matching method */
828be712e4SBarry Smith   redo     = party->redo ? "w3" : "";  /* matching optimization method */
838be712e4SBarry Smith 
848be712e4SBarry Smith   PetscCall(PetscMalloc1(mat->rmap->N, &part_party));
858be712e4SBarry Smith 
868be712e4SBarry Smith   /* redirect output to buffer */
878be712e4SBarry Smith #if defined(PETSC_HAVE_UNISTD_H)
888be712e4SBarry Smith   fd_stdout = dup(1);
898be712e4SBarry Smith   PetscCheck(!pipe(fd_pipe), PETSC_COMM_SELF, PETSC_ERR_SYS, "Could not open pipe");
908be712e4SBarry Smith   close(1);
918be712e4SBarry Smith   dup2(fd_pipe[1], 1);
928be712e4SBarry Smith   PetscCall(PetscMalloc1(SIZE_LOG, &mesg_log));
938be712e4SBarry Smith #endif
948be712e4SBarry Smith 
958be712e4SBarry Smith   /* library call */
968be712e4SBarry Smith   party_lib_times_start();
978be712e4SBarry 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);
988be712e4SBarry Smith 
998be712e4SBarry Smith   party_lib_times_output(1);
1008be712e4SBarry Smith   part_info(n, vertex_w, edge_p, edge, NULL, p, part_party, 1);
1018be712e4SBarry Smith 
1028be712e4SBarry Smith #if defined(PETSC_HAVE_UNISTD_H)
1038be712e4SBarry Smith   PetscCall(PetscFFlush(stdout));
1048be712e4SBarry Smith   count = read(fd_pipe[0], mesg_log, (SIZE_LOG - 1) * sizeof(char));
1058be712e4SBarry Smith   if (count < 0) count = 0;
1068be712e4SBarry Smith   mesg_log[count] = 0;
1078be712e4SBarry Smith   close(1);
1088be712e4SBarry Smith   dup2(fd_stdout, 1);
1098be712e4SBarry Smith   close(fd_stdout);
1108be712e4SBarry Smith   close(fd_pipe[0]);
1118be712e4SBarry Smith   close(fd_pipe[1]);
1128be712e4SBarry Smith   if (party->verbose) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)mat), "%s", mesg_log));
1138be712e4SBarry Smith   PetscCall(PetscFree(mesg_log));
1148be712e4SBarry Smith #endif
1158be712e4SBarry Smith   PetscCheck(!perr, PETSC_COMM_SELF, PETSC_ERR_LIB, "Party failed");
1168be712e4SBarry Smith 
1178be712e4SBarry Smith   PetscCall(PetscMalloc1(mat->rmap->N, &parttab));
1188be712e4SBarry Smith   for (i = 0; i < mat->rmap->N; i++) parttab[i] = part_party[i];
1198be712e4SBarry Smith 
1208be712e4SBarry Smith   /* creation of the index set */
1218be712e4SBarry Smith   nb_locals = mat->rmap->n;
1228be712e4SBarry Smith   locals    = parttab + mat->rmap->rstart;
1238be712e4SBarry Smith 
1248be712e4SBarry Smith   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), nb_locals, locals, PETSC_COPY_VALUES, partitioning));
1258be712e4SBarry Smith 
1268be712e4SBarry Smith   /* clean up */
1278be712e4SBarry Smith   PetscCall(PetscFree(parttab));
1288be712e4SBarry Smith   PetscCall(PetscFree(part_party));
1298be712e4SBarry Smith   PetscCall(MatDestroy(&matSeq));
1308be712e4SBarry Smith   PetscCall(MatDestroy(&matAdj));
1318be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1328be712e4SBarry Smith }
1338be712e4SBarry Smith 
1348be712e4SBarry Smith static PetscErrorCode MatPartitioningView_Party(MatPartitioning part, PetscViewer viewer)
1358be712e4SBarry Smith {
1368be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
1378be712e4SBarry Smith   PetscBool              isascii;
1388be712e4SBarry Smith 
1398be712e4SBarry Smith   PetscFunctionBegin;
1408be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1418be712e4SBarry Smith   if (isascii) {
1428be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "  Global method: %s\n", party->global));
1438be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "  Local method: %s\n", party->local));
1448be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, "  Number of vertices for the coarse graph: %d\n", party->nbvtxcoarsed));
1458be712e4SBarry Smith     if (party->redm) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using matching method for graph reduction\n"));
1468be712e4SBarry Smith     if (party->redo) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using matching optimization\n"));
1478be712e4SBarry Smith     if (party->recursive) PetscCall(PetscViewerASCIIPrintf(viewer, "  Using recursive bipartitioning\n"));
1488be712e4SBarry Smith   }
1498be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1508be712e4SBarry Smith }
1518be712e4SBarry Smith 
1528be712e4SBarry Smith /*@C
1538be712e4SBarry Smith   MatPartitioningPartySetGlobal - Set global method for Party partitioner.
1548be712e4SBarry Smith 
1558be712e4SBarry Smith   Collective
1568be712e4SBarry Smith 
1578be712e4SBarry Smith   Input Parameters:
1588be712e4SBarry Smith + part   - the partitioning context
1598be712e4SBarry Smith - global - a string representing the method
1608be712e4SBarry Smith 
1618be712e4SBarry Smith   Options Database Key:
1628be712e4SBarry Smith . -mat_partitioning_party_global <method> - the global method
1638be712e4SBarry Smith 
1648be712e4SBarry Smith   Level: advanced
1658be712e4SBarry Smith 
1668be712e4SBarry Smith   Note:
1678be712e4SBarry Smith   The method may be one of `MP_PARTY_OPT`, `MP_PARTY_LIN`, `MP_PARTY_SCA`,
1688be712e4SBarry Smith   `MP_PARTY_RAN`, `MP_PARTY_GBF`, `MP_PARTY_GCF`, `MP_PARTY_BUB` or `MP_PARTY_DEF`, or
1698be712e4SBarry Smith   alternatively a string describing the method. Two or more methods can be
1708be712e4SBarry Smith   combined like "gbf,gcf". Check the Party Library Users Manual for details.
1718be712e4SBarry Smith 
1728be712e4SBarry Smith .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetLocal()`
1738be712e4SBarry Smith @*/
1748be712e4SBarry Smith PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part, const char *global)
1758be712e4SBarry Smith {
1768be712e4SBarry Smith   PetscFunctionBegin;
1778be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
1788be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPartySetGlobal_C", (MatPartitioning, const char *), (part, global));
1798be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1808be712e4SBarry Smith }
1818be712e4SBarry Smith 
1828be712e4SBarry Smith static PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part, const char *global)
1838be712e4SBarry Smith {
1848be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
1858be712e4SBarry Smith 
1868be712e4SBarry Smith   PetscFunctionBegin;
1878be712e4SBarry Smith   PetscCall(PetscStrncpy(party->global, global, 15));
1888be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
1898be712e4SBarry Smith }
1908be712e4SBarry Smith 
1918be712e4SBarry Smith /*@C
1928be712e4SBarry Smith   MatPartitioningPartySetLocal - Set local method used by the Party partitioner.
1938be712e4SBarry Smith 
1948be712e4SBarry Smith   Collective
1958be712e4SBarry Smith 
1968be712e4SBarry Smith   Input Parameters:
1978be712e4SBarry Smith + part  - the partitioning context
1988be712e4SBarry Smith - local - a string representing the method
1998be712e4SBarry Smith 
2008be712e4SBarry Smith   Options Database Key:
2018be712e4SBarry Smith . -mat_partitioning_party_local <method> - the local method
2028be712e4SBarry Smith 
2038be712e4SBarry Smith   Level: advanced
2048be712e4SBarry Smith 
2058be712e4SBarry Smith   Note:
2068be712e4SBarry Smith   The method may be one of `MP_PARTY_HELPFUL_SETS`, `MP_PARTY_KERNIGHAN_LIN`, or
2078be712e4SBarry Smith   `MP_PARTY_NONE`. Check the Party Library Users Manual for details.
2088be712e4SBarry Smith 
2098be712e4SBarry Smith .seealso: `MATPARTITIONINGPARTY`, `MatPartitioningPartySetGlobal()`
2108be712e4SBarry Smith @*/
2118be712e4SBarry Smith PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part, const char *local)
2128be712e4SBarry Smith {
2138be712e4SBarry Smith   PetscFunctionBegin;
2148be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
2158be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPartySetLocal_C", (MatPartitioning, const char *), (part, local));
2168be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2178be712e4SBarry Smith }
2188be712e4SBarry Smith 
2198be712e4SBarry Smith static PetscErrorCode MatPartitioningPartySetLocal_Party(MatPartitioning part, const char *local)
2208be712e4SBarry Smith 
2218be712e4SBarry Smith {
2228be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
2238be712e4SBarry Smith 
2248be712e4SBarry Smith   PetscFunctionBegin;
2258be712e4SBarry Smith   PetscCall(PetscStrncpy(party->local, local, 15));
2268be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2278be712e4SBarry Smith }
2288be712e4SBarry Smith 
2298be712e4SBarry Smith /*@
2308be712e4SBarry Smith   MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
2318be712e4SBarry Smith   Party partitioner.
2328be712e4SBarry Smith 
2338be712e4SBarry Smith   Collective
2348be712e4SBarry Smith 
2358be712e4SBarry Smith   Input Parameters:
2368be712e4SBarry Smith + part  - the partitioning context
2378be712e4SBarry Smith - level - the coarse level in range [0.0,1.0]
2388be712e4SBarry Smith 
2398be712e4SBarry Smith   Options Database Key:
2408be712e4SBarry Smith . -mat_partitioning_party_coarse <l> - Coarse level
2418be712e4SBarry Smith 
2428be712e4SBarry Smith   Level: advanced
2438be712e4SBarry Smith 
2448be712e4SBarry Smith .seealso: `MATPARTITIONINGPARTY`
2458be712e4SBarry Smith @*/
2468be712e4SBarry Smith PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part, PetscReal level)
2478be712e4SBarry Smith {
2488be712e4SBarry Smith   PetscFunctionBegin;
2498be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
2508be712e4SBarry Smith   PetscValidLogicalCollectiveReal(part, level, 2);
2518be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPartySetCoarseLevel_C", (MatPartitioning, PetscReal), (part, level));
2528be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2538be712e4SBarry Smith }
2548be712e4SBarry Smith 
2558be712e4SBarry Smith static PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part, PetscReal level)
2568be712e4SBarry Smith {
2578be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
2588be712e4SBarry Smith 
2598be712e4SBarry Smith   PetscFunctionBegin;
2608be712e4SBarry 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]");
2618be712e4SBarry Smith   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
2628be712e4SBarry Smith   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
2638be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2648be712e4SBarry Smith }
2658be712e4SBarry Smith 
2668be712e4SBarry Smith /*@
2678be712e4SBarry Smith   MatPartitioningPartySetMatchOptimization - Activate matching optimization for
2688be712e4SBarry Smith   graph reduction.
2698be712e4SBarry Smith 
2708be712e4SBarry Smith   Collective
2718be712e4SBarry Smith 
2728be712e4SBarry Smith   Input Parameters:
2738be712e4SBarry Smith + part - the partitioning context
2748be712e4SBarry Smith - opt  - boolean flag
2758be712e4SBarry Smith 
2768be712e4SBarry Smith   Options Database Key:
2778be712e4SBarry Smith . -mat_partitioning_party_match_optimization - Matching optimization on/off
2788be712e4SBarry Smith 
2798be712e4SBarry Smith   Level: advanced
2808be712e4SBarry Smith 
2818be712e4SBarry Smith .seealso: `MATPARTITIONINGPARTY`
2828be712e4SBarry Smith @*/
2838be712e4SBarry Smith PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part, PetscBool opt)
2848be712e4SBarry Smith {
2858be712e4SBarry Smith   PetscFunctionBegin;
2868be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
2878be712e4SBarry Smith   PetscValidLogicalCollectiveBool(part, opt, 2);
2888be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPartySetMatchOptimization_C", (MatPartitioning, PetscBool), (part, opt));
2898be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2908be712e4SBarry Smith }
2918be712e4SBarry Smith 
2928be712e4SBarry Smith static PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part, PetscBool opt)
2938be712e4SBarry Smith {
2948be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
2958be712e4SBarry Smith 
2968be712e4SBarry Smith   PetscFunctionBegin;
2978be712e4SBarry Smith   party->redo = opt;
2988be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
2998be712e4SBarry Smith }
3008be712e4SBarry Smith 
3018be712e4SBarry Smith /*@
3028be712e4SBarry Smith   MatPartitioningPartySetBipart - Activate or deactivate recursive bisection in the Party partitioner
3038be712e4SBarry Smith 
3048be712e4SBarry Smith   Collective
3058be712e4SBarry Smith 
3068be712e4SBarry Smith   Input Parameters:
3078be712e4SBarry Smith + part - the partitioning context
3088be712e4SBarry Smith - bp   - boolean flag
3098be712e4SBarry Smith 
3108be712e4SBarry Smith   Options Database Key:
3118be712e4SBarry Smith . -mat_partitioning_party_bipart - Bipartitioning option on/off
3128be712e4SBarry Smith 
3138be712e4SBarry Smith   Level: advanced
3148be712e4SBarry Smith 
3158be712e4SBarry Smith .seealso: `MATPARTITIONINGPARTY`
3168be712e4SBarry Smith @*/
3178be712e4SBarry Smith PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part, PetscBool bp)
3188be712e4SBarry Smith {
3198be712e4SBarry Smith   PetscFunctionBegin;
3208be712e4SBarry Smith   PetscValidHeaderSpecific(part, MAT_PARTITIONING_CLASSID, 1);
3218be712e4SBarry Smith   PetscValidLogicalCollectiveBool(part, bp, 2);
3228be712e4SBarry Smith   PetscTryMethod(part, "MatPartitioningPartySetBipart_C", (MatPartitioning, PetscBool), (part, bp));
3238be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3248be712e4SBarry Smith }
3258be712e4SBarry Smith 
3268be712e4SBarry Smith static PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part, PetscBool bp)
3278be712e4SBarry Smith {
3288be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
3298be712e4SBarry Smith 
3308be712e4SBarry Smith   PetscFunctionBegin;
3318be712e4SBarry Smith   party->recursive = bp;
3328be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3338be712e4SBarry Smith }
3348be712e4SBarry Smith 
3358be712e4SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_Party(MatPartitioning part, PetscOptionItems *PetscOptionsObject)
3368be712e4SBarry Smith {
3378be712e4SBarry Smith   PetscBool              flag;
3388be712e4SBarry Smith   char                   value[256];
3398be712e4SBarry Smith   PetscReal              r;
3408be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
3418be712e4SBarry Smith 
3428be712e4SBarry Smith   PetscFunctionBegin;
3438be712e4SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Set Party partitioning options");
3448be712e4SBarry Smith   PetscCall(PetscOptionsString("-mat_partitioning_party_global", "Global method", "MatPartitioningPartySetGlobal", party->global, value, sizeof(value), &flag));
3458be712e4SBarry Smith   if (flag) PetscCall(MatPartitioningPartySetGlobal(part, value));
3468be712e4SBarry Smith   PetscCall(PetscOptionsString("-mat_partitioning_party_local", "Local method", "MatPartitioningPartySetLocal", party->local, value, sizeof(value), &flag));
3478be712e4SBarry Smith   if (flag) PetscCall(MatPartitioningPartySetLocal(part, value));
3488be712e4SBarry Smith   PetscCall(PetscOptionsReal("-mat_partitioning_party_coarse", "Coarse level", "MatPartitioningPartySetCoarseLevel", 0.0, &r, &flag));
3498be712e4SBarry Smith   if (flag) PetscCall(MatPartitioningPartySetCoarseLevel(part, r));
3508be712e4SBarry Smith   PetscCall(PetscOptionsBool("-mat_partitioning_party_match_optimization", "Matching optimization on/off", "MatPartitioningPartySetMatchOptimization", party->redo, &party->redo, NULL));
3518be712e4SBarry Smith   PetscCall(PetscOptionsBool("-mat_partitioning_party_bipart", "Bipartitioning on/off", "MatPartitioningPartySetBipart", party->recursive, &party->recursive, NULL));
3528be712e4SBarry Smith   PetscCall(PetscOptionsBool("-mat_partitioning_party_verbose", "Show library output", "", party->verbose, &party->verbose, NULL));
3538be712e4SBarry Smith   PetscOptionsHeadEnd();
3548be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3558be712e4SBarry Smith }
3568be712e4SBarry Smith 
3578be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
3588be712e4SBarry Smith {
3598be712e4SBarry Smith   MatPartitioning_Party *party = (MatPartitioning_Party *)part->data;
3608be712e4SBarry Smith 
3618be712e4SBarry Smith   PetscFunctionBegin;
3628be712e4SBarry Smith   PetscCall(PetscFree(party));
3638be712e4SBarry Smith   /* clear composed functions */
3648be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", NULL));
3658be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", NULL));
3668be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", NULL));
3678be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", NULL));
3688be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", NULL));
3698be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
3708be712e4SBarry Smith }
3718be712e4SBarry Smith 
3728be712e4SBarry Smith /*MC
373*613ce9feSSatish Balay    MATPARTITIONINGPARTY - Creates a partitioning context via the external package Party <
374*613ce9feSSatish Balay    http://wwwcs.upb.de/fachbereich/AG/monien/RESEARCH/PART/party.htm>.
3758be712e4SBarry Smith 
3768be712e4SBarry Smith    Level: beginner
3778be712e4SBarry Smith 
378*613ce9feSSatish Balay    Note:
3798be712e4SBarry Smith    Does not support the `MatPartitioningSetUseEdgeWeights()` option
3808be712e4SBarry Smith 
3818be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MatPartitioningPartySetGlobal()`, `MatPartitioningPartySetLocal()`,
3828be712e4SBarry Smith           `MatPartitioningPartySetCoarseLevel()`, `MatPartitioningPartySetMatchOptimization()`, `MatPartitioningPartySetBipart()`
3838be712e4SBarry Smith M*/
3848be712e4SBarry Smith 
3858be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning part)
3868be712e4SBarry Smith {
3878be712e4SBarry Smith   MatPartitioning_Party *party;
3888be712e4SBarry Smith 
3898be712e4SBarry Smith   PetscFunctionBegin;
3908be712e4SBarry Smith   PetscCall(PetscNew(&party));
3918be712e4SBarry Smith   part->data = (void *)party;
3928be712e4SBarry Smith 
3938be712e4SBarry Smith   PetscCall(PetscStrncpy(party->global, "gcf,gbf", sizeof(party->global)));
3948be712e4SBarry Smith   PetscCall(PetscStrncpy(party->local, "kl", sizeof(party->local)));
3958be712e4SBarry Smith 
3968be712e4SBarry Smith   party->redm         = PETSC_TRUE;
3978be712e4SBarry Smith   party->redo         = PETSC_TRUE;
3988be712e4SBarry Smith   party->recursive    = PETSC_TRUE;
3998be712e4SBarry Smith   party->verbose      = PETSC_FALSE;
4008be712e4SBarry Smith   party->nbvtxcoarsed = 200;
4018be712e4SBarry Smith 
4028be712e4SBarry Smith   part->ops->apply          = MatPartitioningApply_Party;
4038be712e4SBarry Smith   part->ops->view           = MatPartitioningView_Party;
4048be712e4SBarry Smith   part->ops->destroy        = MatPartitioningDestroy_Party;
4058be712e4SBarry Smith   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;
4068be712e4SBarry Smith 
4078be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetGlobal_C", MatPartitioningPartySetGlobal_Party));
4088be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetLocal_C", MatPartitioningPartySetLocal_Party));
4098be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetCoarseLevel_C", MatPartitioningPartySetCoarseLevel_Party));
4108be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetMatchOptimization_C", MatPartitioningPartySetMatchOptimization_Party));
4118be712e4SBarry Smith   PetscCall(PetscObjectComposeFunction((PetscObject)part, "MatPartitioningPartySetBipart_C", MatPartitioningPartySetBipart_Party));
4128be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
4138be712e4SBarry Smith }
414