xref: /petsc/src/dm/partitioner/impls/parmetis/partparmetis.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
2abe9303eSLisandro Dalcin 
3abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PARMETIS)
4abe9303eSLisandro Dalcin   #include <parmetis.h>
5abe9303eSLisandro Dalcin #endif
6abe9303eSLisandro Dalcin 
7abe9303eSLisandro Dalcin PetscBool  ParMetisPartitionerCite       = PETSC_FALSE;
89371c9d4SSatish Balay const char ParMetisPartitionerCitation[] = "@article{KarypisKumar98,\n"
9abe9303eSLisandro Dalcin                                            "  author  = {George Karypis and Vipin Kumar},\n"
10abe9303eSLisandro Dalcin                                            "  title   = {A Parallel Algorithm for Multilevel Graph Partitioning and Sparse Matrix Ordering},\n"
11abe9303eSLisandro Dalcin                                            "  journal = {Journal of Parallel and Distributed Computing},\n"
12abe9303eSLisandro Dalcin                                            "  volume  = {48},\n"
13abe9303eSLisandro Dalcin                                            "  pages   = {71--85},\n"
14abe9303eSLisandro Dalcin                                            "  year    = {1998}\n"
15abe9303eSLisandro Dalcin                                            "  doi     = {https://doi.org/10.1006/jpdc.1997.1403}\n"
16abe9303eSLisandro Dalcin                                            "}\n";
17abe9303eSLisandro Dalcin 
18abe9303eSLisandro Dalcin typedef struct {
19abe9303eSLisandro Dalcin   MPI_Comm  pcomm;
20abe9303eSLisandro Dalcin   PetscInt  ptype;
21abe9303eSLisandro Dalcin   PetscReal imbalanceRatio;
22abe9303eSLisandro Dalcin   PetscInt  debugFlag;
23abe9303eSLisandro Dalcin   PetscInt  randomSeed;
24abe9303eSLisandro Dalcin } PetscPartitioner_ParMetis;
25abe9303eSLisandro Dalcin 
26abe9303eSLisandro Dalcin static const char *ptypes[] = {"kway", "rb"};
27abe9303eSLisandro Dalcin 
PetscPartitionerDestroy_ParMetis(PetscPartitioner part)28d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part)
29d71ae5a4SJacob Faibussowitsch {
30abe9303eSLisandro Dalcin   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *)part->data;
31abe9303eSLisandro Dalcin 
32abe9303eSLisandro Dalcin   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&p->pcomm));
349566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36abe9303eSLisandro Dalcin }
37abe9303eSLisandro Dalcin 
PetscPartitionerView_ParMetis_ASCII(PetscPartitioner part,PetscViewer viewer)38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_ParMetis_ASCII(PetscPartitioner part, PetscViewer viewer)
39d71ae5a4SJacob Faibussowitsch {
40abe9303eSLisandro Dalcin   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *)part->data;
41abe9303eSLisandro Dalcin 
42abe9303eSLisandro Dalcin   PetscFunctionBegin;
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
449566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "ParMetis type: %s\n", ptypes[p->ptype]));
459566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "load imbalance ratio %g\n", (double)p->imbalanceRatio));
4663a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "debug flag %" PetscInt_FMT "\n", p->debugFlag));
4763a3b9bcSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "random seed %" PetscInt_FMT "\n", p->randomSeed));
489566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50abe9303eSLisandro Dalcin }
51abe9303eSLisandro Dalcin 
PetscPartitionerView_ParMetis(PetscPartitioner part,PetscViewer viewer)52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer)
53d71ae5a4SJacob Faibussowitsch {
54*9f196a02SMartin Diehl   PetscBool isascii;
55abe9303eSLisandro Dalcin 
56abe9303eSLisandro Dalcin   PetscFunctionBegin;
57abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
58abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
59*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
60*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscPartitionerView_ParMetis_ASCII(part, viewer));
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
62abe9303eSLisandro Dalcin }
63abe9303eSLisandro Dalcin 
PetscPartitionerSetFromOptions_ParMetis(PetscPartitioner part,PetscOptionItems PetscOptionsObject)64ce78bad3SBarry Smith static PetscErrorCode PetscPartitionerSetFromOptions_ParMetis(PetscPartitioner part, PetscOptionItems PetscOptionsObject)
65d71ae5a4SJacob Faibussowitsch {
66abe9303eSLisandro Dalcin   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *)part->data;
67abe9303eSLisandro Dalcin 
68abe9303eSLisandro Dalcin   PetscFunctionBegin;
69d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner ParMetis Options");
709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-petscpartitioner_parmetis_type", "Partitioning method", "", ptypes, 2, ptypes[p->ptype], &p->ptype, NULL));
719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-petscpartitioner_parmetis_imbalance_ratio", "Load imbalance ratio limit", "", p->imbalanceRatio, &p->imbalanceRatio, NULL));
729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-petscpartitioner_parmetis_debug", "Debugging flag", "", p->debugFlag, &p->debugFlag, NULL));
739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-petscpartitioner_parmetis_seed", "Random seed", "", p->randomSeed, &p->randomSeed, NULL));
74d0609cedSBarry Smith   PetscOptionsHeadEnd();
753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
76abe9303eSLisandro Dalcin }
77abe9303eSLisandro Dalcin 
PetscPartitionerPartition_ParMetis(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * partition)7821c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
79d71ae5a4SJacob Faibussowitsch {
80abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PARMETIS)
81abe9303eSLisandro Dalcin   PetscPartitioner_ParMetis *pm = (PetscPartitioner_ParMetis *)part->data;
82abe9303eSLisandro Dalcin   MPI_Comm                   comm;
83abe9303eSLisandro Dalcin   PetscInt                   nvtxs = numVertices;     /* The number of vertices in full graph */
84abe9303eSLisandro Dalcin   PetscInt                  *vtxdist;                 /* Distribution of vertices across processes */
85abe9303eSLisandro Dalcin   PetscInt                  *xadj        = start;     /* Start of edge list for each vertex */
86abe9303eSLisandro Dalcin   PetscInt                  *adjncy      = adjacency; /* Edge lists for all vertices */
87abe9303eSLisandro Dalcin   PetscInt                  *vwgt        = NULL;      /* Vertex weights */
88abe9303eSLisandro Dalcin   PetscInt                  *adjwgt      = NULL;      /* Edge weights */
89abe9303eSLisandro Dalcin   PetscInt                   wgtflag     = 0;         /* Indicates which weights are present */
90abe9303eSLisandro Dalcin   PetscInt                   numflag     = 0;         /* Indicates initial offset (0 or 1) */
91abe9303eSLisandro Dalcin   PetscInt                   ncon        = 1;         /* The number of weights per vertex */
92abe9303eSLisandro Dalcin   PetscInt                   metis_ptype = pm->ptype; /* kway or recursive bisection */
93abe9303eSLisandro Dalcin   real_t                    *tpwgts;                  /* The fraction of vertex weights assigned to each partition */
94abe9303eSLisandro Dalcin   real_t                    *ubvec;                   /* The balance intolerance for vertex weights */
95abe9303eSLisandro Dalcin   PetscInt                   options[64];             /* Options */
96abe9303eSLisandro Dalcin   PetscInt                   v, i, *assignment, *points;
97abe9303eSLisandro Dalcin   PetscMPIInt                p, size, rank;
98abe9303eSLisandro Dalcin   PetscBool                  hasempty = PETSC_FALSE;
99abe9303eSLisandro Dalcin 
100abe9303eSLisandro Dalcin   PetscFunctionBegin;
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
1029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
104abe9303eSLisandro Dalcin   /* Calculate vertex distribution */
1059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc4(size + 1, &vtxdist, nparts * ncon, &tpwgts, ncon, &ubvec, nvtxs, &assignment));
106abe9303eSLisandro Dalcin   vtxdist[0] = 0;
1079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm));
108abe9303eSLisandro Dalcin   for (p = 2; p <= size; ++p) {
109abe9303eSLisandro Dalcin     hasempty = (PetscBool)(hasempty || !vtxdist[p - 1] || !vtxdist[p]);
110abe9303eSLisandro Dalcin     vtxdist[p] += vtxdist[p - 1];
111abe9303eSLisandro Dalcin   }
112abe9303eSLisandro Dalcin   /* null graph */
113abe9303eSLisandro Dalcin   if (vtxdist[size] == 0) {
1149566063dSJacob Faibussowitsch     PetscCall(PetscFree4(vtxdist, tpwgts, ubvec, assignment));
1159566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition));
1163ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
117abe9303eSLisandro Dalcin   }
118abe9303eSLisandro Dalcin   /* Calculate partition weights */
119abe9303eSLisandro Dalcin   if (targetSection) {
120abe9303eSLisandro Dalcin     PetscInt p;
121abe9303eSLisandro Dalcin     real_t   sumt = 0.0;
122abe9303eSLisandro Dalcin 
123abe9303eSLisandro Dalcin     for (p = 0; p < nparts; ++p) {
124abe9303eSLisandro Dalcin       PetscInt tpd;
125abe9303eSLisandro Dalcin 
1269566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(targetSection, p, &tpd));
127abe9303eSLisandro Dalcin       sumt += tpd;
128abe9303eSLisandro Dalcin       tpwgts[p] = tpd;
129abe9303eSLisandro Dalcin     }
130abe9303eSLisandro Dalcin     if (sumt) { /* METIS/ParMETIS do not like exactly zero weight */
131abe9303eSLisandro Dalcin       for (p = 0, sumt = 0.0; p < nparts; ++p) {
1326497c311SBarry Smith         tpwgts[p] = (real_t)PetscMax(tpwgts[p], PETSC_SMALL);
133abe9303eSLisandro Dalcin         sumt += tpwgts[p];
134abe9303eSLisandro Dalcin       }
135abe9303eSLisandro Dalcin       for (p = 0; p < nparts; ++p) tpwgts[p] /= sumt;
136abe9303eSLisandro Dalcin       for (p = 0, sumt = 0.0; p < nparts - 1; ++p) sumt += tpwgts[p];
1376497c311SBarry Smith       tpwgts[nparts - 1] = (real_t)(1. - sumt);
138abe9303eSLisandro Dalcin     }
139abe9303eSLisandro Dalcin   } else {
1406497c311SBarry Smith     for (p = 0; p < nparts; ++p) tpwgts[p] = (real_t)(1.0 / nparts);
141abe9303eSLisandro Dalcin   }
1426497c311SBarry Smith   ubvec[0] = (real_t)pm->imbalanceRatio;
143abe9303eSLisandro Dalcin 
144abe9303eSLisandro Dalcin   /* Weight cells */
145abe9303eSLisandro Dalcin   if (vertSection) {
1469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nvtxs, &vwgt));
14748a46eb9SPierre Jolivet     for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionGetDof(vertSection, v, &vwgt[v]));
148abe9303eSLisandro Dalcin     wgtflag |= 2; /* have weights on graph vertices */
149abe9303eSLisandro Dalcin   }
15021c92275SMatthew G. Knepley   // Weight edges
15121c92275SMatthew G. Knepley   if (edgeSection) {
15221c92275SMatthew G. Knepley     PetscCall(PetscMalloc1(xadj[nvtxs], &adjwgt));
15321c92275SMatthew G. Knepley     for (PetscInt e = 0; e < xadj[nvtxs]; ++e) PetscCall(PetscSectionGetDof(edgeSection, e, &adjwgt[e]));
15421c92275SMatthew G. Knepley     wgtflag |= 1; /* have weights on graph edges */
15521c92275SMatthew G. Knepley   }
156abe9303eSLisandro Dalcin 
157fbccb6d4SPierre Jolivet   for (p = 0; !vtxdist[p + 1] && p < size; ++p);
158abe9303eSLisandro Dalcin   if (vtxdist[p + 1] == vtxdist[size]) {
159abe9303eSLisandro Dalcin     if (rank == p) {
1602da392ccSBarry Smith       int err;
1612da392ccSBarry Smith       err                          = METIS_SetDefaultOptions(options); /* initialize all defaults */
162abe9303eSLisandro Dalcin       options[METIS_OPTION_DBGLVL] = pm->debugFlag;
163abe9303eSLisandro Dalcin       options[METIS_OPTION_SEED]   = pm->randomSeed;
16408401ef6SPierre Jolivet       PetscCheck(err == METIS_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_SetDefaultOptions()");
165abe9303eSLisandro Dalcin       if (metis_ptype == 1) {
166792fecdfSBarry Smith         PetscStackPushExternal("METIS_PartGraphRecursive");
1672da392ccSBarry Smith         err = METIS_PartGraphRecursive(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, options, &part->edgeCut, assignment);
168abe9303eSLisandro Dalcin         PetscStackPop;
16908401ef6SPierre Jolivet         PetscCheck(err == METIS_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphRecursive()");
170abe9303eSLisandro Dalcin       } else {
171abe9303eSLisandro Dalcin         /*
172abe9303eSLisandro Dalcin          It would be nice to activate the two options below, but they would need some actual testing.
173abe9303eSLisandro Dalcin          - Turning on these options may exercise path of the METIS code that have bugs and may break production runs.
174abe9303eSLisandro Dalcin          - If CONTIG is set to 1, METIS will exit with error if the graph is disconnected, despite the manual saying the option is ignored in such case.
175abe9303eSLisandro Dalcin         */
176abe9303eSLisandro Dalcin         /* options[METIS_OPTION_CONTIG]  = 1; */ /* try to produce partitions that are contiguous */
177abe9303eSLisandro Dalcin         /* options[METIS_OPTION_MINCONN] = 1; */ /* minimize the maximum degree of the subdomain graph */
178792fecdfSBarry Smith         PetscStackPushExternal("METIS_PartGraphKway");
1792da392ccSBarry Smith         err = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, options, &part->edgeCut, assignment);
180abe9303eSLisandro Dalcin         PetscStackPop;
18108401ef6SPierre Jolivet         PetscCheck(err == METIS_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
182abe9303eSLisandro Dalcin       }
183abe9303eSLisandro Dalcin     }
184abe9303eSLisandro Dalcin   } else {
185abe9303eSLisandro Dalcin     MPI_Comm pcomm = pm->pcomm;
186abe9303eSLisandro Dalcin 
187abe9303eSLisandro Dalcin     options[0] = 1; /*use options */
188abe9303eSLisandro Dalcin     options[1] = pm->debugFlag;
189abe9303eSLisandro Dalcin     options[2] = (pm->randomSeed == -1) ? 15 : pm->randomSeed; /* default is GLOBAL_SEED=15 from `libparmetis/defs.h` */
190abe9303eSLisandro Dalcin 
191abe9303eSLisandro Dalcin     if (hasempty) { /* parmetis does not support empty graphs on some of the processes */
192abe9303eSLisandro Dalcin       PetscInt cnt;
193abe9303eSLisandro Dalcin 
1949566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_split(pm->pcomm, !!nvtxs, rank, &pcomm));
195abe9303eSLisandro Dalcin       for (p = 0, cnt = 0; p < size; p++) {
196abe9303eSLisandro Dalcin         if (vtxdist[p + 1] != vtxdist[p]) {
197abe9303eSLisandro Dalcin           vtxdist[cnt + 1] = vtxdist[p + 1];
198abe9303eSLisandro Dalcin           cnt++;
199abe9303eSLisandro Dalcin         }
200abe9303eSLisandro Dalcin       }
201abe9303eSLisandro Dalcin     }
202abe9303eSLisandro Dalcin     if (nvtxs) {
2032da392ccSBarry Smith       int err;
204792fecdfSBarry Smith       PetscStackPushExternal("ParMETIS_V3_PartKway");
2052da392ccSBarry Smith       err = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &part->edgeCut, assignment, &pcomm);
206abe9303eSLisandro Dalcin       PetscStackPop;
20708401ef6SPierre Jolivet       PetscCheck(err == METIS_OK, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error %d in ParMETIS_V3_PartKway()", err);
208abe9303eSLisandro Dalcin     }
20948a46eb9SPierre Jolivet     if (hasempty) PetscCallMPI(MPI_Comm_free(&pcomm));
210abe9303eSLisandro Dalcin   }
211abe9303eSLisandro Dalcin 
212abe9303eSLisandro Dalcin   /* Convert to PetscSection+IS */
2139566063dSJacob Faibussowitsch   for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1));
2149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nvtxs, &points));
215abe9303eSLisandro Dalcin   for (p = 0, i = 0; p < nparts; ++p) {
216abe9303eSLisandro Dalcin     for (v = 0; v < nvtxs; ++v) {
217abe9303eSLisandro Dalcin       if (assignment[v] == p) points[i++] = v;
218abe9303eSLisandro Dalcin     }
219abe9303eSLisandro Dalcin   }
22063a3b9bcSJacob Faibussowitsch   PetscCheck(i == nvtxs, comm, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs);
2219566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition));
2229566063dSJacob Faibussowitsch   PetscCall(PetscFree4(vtxdist, tpwgts, ubvec, assignment));
2239566063dSJacob Faibussowitsch   PetscCall(PetscFree(vwgt));
22421c92275SMatthew G. Knepley   PetscCall(PetscFree(adjwgt));
2253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
226abe9303eSLisandro Dalcin #else
227abe9303eSLisandro Dalcin   SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis.");
228abe9303eSLisandro Dalcin #endif
229abe9303eSLisandro Dalcin }
230abe9303eSLisandro Dalcin 
PetscPartitionerInitialize_ParMetis(PetscPartitioner part)231d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part)
232d71ae5a4SJacob Faibussowitsch {
233abe9303eSLisandro Dalcin   PetscFunctionBegin;
234abe9303eSLisandro Dalcin   part->noGraph             = PETSC_FALSE;
235abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_ParMetis;
236abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_ParMetis;
237abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_ParMetis;
238abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_ParMetis;
2393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
240abe9303eSLisandro Dalcin }
241abe9303eSLisandro Dalcin 
242abe9303eSLisandro Dalcin /*MC
243abe9303eSLisandro Dalcin   PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMETIS library
244abe9303eSLisandro Dalcin 
245abe9303eSLisandro Dalcin   Level: intermediate
246abe9303eSLisandro Dalcin 
247abe9303eSLisandro Dalcin   Options Database Keys:
248abe9303eSLisandro Dalcin +  -petscpartitioner_parmetis_type <string> - ParMETIS partitioning type. Either "kway" or "rb" (recursive bisection)
249abe9303eSLisandro Dalcin .  -petscpartitioner_parmetis_imbalance_ratio <value> - Load imbalance ratio limit
250abe9303eSLisandro Dalcin .  -petscpartitioner_parmetis_debug <int> - Debugging flag passed to ParMETIS/METIS routines
251abe9303eSLisandro Dalcin -  -petscpartitioner_parmetis_seed <int> - Random seed
252abe9303eSLisandro Dalcin 
253abe9303eSLisandro Dalcin   Notes: when the graph is on a single process, this partitioner actually calls METIS and not ParMETIS
254abe9303eSLisandro Dalcin 
255db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
256abe9303eSLisandro Dalcin M*/
257abe9303eSLisandro Dalcin 
PetscPartitionerCreate_ParMetis(PetscPartitioner part)258d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part)
259d71ae5a4SJacob Faibussowitsch {
260abe9303eSLisandro Dalcin   PetscPartitioner_ParMetis *p;
261abe9303eSLisandro Dalcin 
262abe9303eSLisandro Dalcin   PetscFunctionBegin;
263abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
2644dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
265abe9303eSLisandro Dalcin   part->data = p;
266abe9303eSLisandro Dalcin 
2679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)part), &p->pcomm));
268abe9303eSLisandro Dalcin   p->ptype          = 0;
269abe9303eSLisandro Dalcin   p->imbalanceRatio = 1.05;
270abe9303eSLisandro Dalcin   p->debugFlag      = 0;
271abe9303eSLisandro Dalcin   p->randomSeed     = -1; /* defaults to GLOBAL_SEED=15 from `libparmetis/defs.h` */
272abe9303eSLisandro Dalcin 
2739566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_ParMetis(part));
2749566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionerCite));
2753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
276abe9303eSLisandro Dalcin }
277