xref: /petsc/src/dm/partitioner/impls/simple/partsimple.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
163a3b9bcSJacob Faibussowitsch 
2abe9303eSLisandro Dalcin #include <petscvec.h>
3abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
4abe9303eSLisandro Dalcin 
5abe9303eSLisandro Dalcin typedef struct {
60c569c6eSMark   PetscBool useGrid;        /* Flag to use a grid layout */
70c569c6eSMark   PetscInt  gridDim;        /* The grid dimension */
80c569c6eSMark   PetscInt  nodeGrid[3];    /* Dimension of node grid */
90c569c6eSMark   PetscInt  processGrid[3]; /* Dimension of local process grid on each node */
10abe9303eSLisandro Dalcin } PetscPartitioner_Simple;
11abe9303eSLisandro Dalcin 
12d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part)
13d71ae5a4SJacob Faibussowitsch {
14abe9303eSLisandro Dalcin   PetscFunctionBegin;
159566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
16*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17abe9303eSLisandro Dalcin }
18abe9303eSLisandro Dalcin 
19d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Simple_ASCII(PetscPartitioner part, PetscViewer viewer)
20d71ae5a4SJacob Faibussowitsch {
21abe9303eSLisandro Dalcin   PetscFunctionBegin;
22*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23abe9303eSLisandro Dalcin }
24abe9303eSLisandro Dalcin 
25d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer)
26d71ae5a4SJacob Faibussowitsch {
27abe9303eSLisandro Dalcin   PetscBool iascii;
28abe9303eSLisandro Dalcin 
29abe9303eSLisandro Dalcin   PetscFunctionBegin;
30abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
31abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
339566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Simple_ASCII(part, viewer));
34*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35abe9303eSLisandro Dalcin }
36abe9303eSLisandro Dalcin 
37d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerSetFromOptions_Simple(PetscPartitioner part, PetscOptionItems *PetscOptionsObject)
38d71ae5a4SJacob Faibussowitsch {
390c569c6eSMark   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *)part->data;
400c569c6eSMark   PetscInt                 num, i;
410c569c6eSMark   PetscBool                flg;
420c569c6eSMark 
430c569c6eSMark   PetscFunctionBegin;
440c569c6eSMark   for (i = 0; i < 3; ++i) p->processGrid[i] = p->nodeGrid[i] = 1;
45d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner Simple Options");
460c569c6eSMark   num = 3;
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_node_grid", "Number of nodes in each dimension", "", p->nodeGrid, &num, &flg));
489371c9d4SSatish Balay   if (flg) {
499371c9d4SSatish Balay     p->useGrid = PETSC_TRUE;
509371c9d4SSatish Balay     p->gridDim = num;
519371c9d4SSatish Balay   }
520c569c6eSMark   num = 3;
539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_process_grid", "Number of local processes in each dimension for a given node", "", p->processGrid, &num, &flg));
540c569c6eSMark   if (flg) {
550c569c6eSMark     p->useGrid = PETSC_TRUE;
560c569c6eSMark     if (p->gridDim < 0) p->gridDim = num;
5763a3b9bcSJacob Faibussowitsch     else PetscCheck(p->gridDim == num, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_INCOMP, "Process grid dimension %" PetscInt_FMT " != %" PetscInt_FMT " node grid dimension", num, p->gridDim);
580c569c6eSMark   }
59d0609cedSBarry Smith   PetscOptionsHeadEnd();
60*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
610c569c6eSMark }
620c569c6eSMark 
63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerPartition_Simple_Grid(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
64d71ae5a4SJacob Faibussowitsch {
650c569c6eSMark   PetscPartitioner_Simple *p     = (PetscPartitioner_Simple *)part->data;
660c569c6eSMark   const PetscInt          *nodes = p->nodeGrid;
670c569c6eSMark   const PetscInt          *procs = p->processGrid;
680c569c6eSMark   PetscInt                *cellproc, *offsets, cells[3] = {1, 1, 1}, pcells[3] = {1, 1, 1};
69bd026e97SJed Brown   PetscInt                 Np = 1, Nr, np, nk, nj, ni, pk, pj, pi, ck, cj, ci, i;
700c569c6eSMark   MPI_Comm                 comm;
710c569c6eSMark   PetscMPIInt              size;
720c569c6eSMark 
730c569c6eSMark   PetscFunctionBegin;
749566063dSJacob Faibussowitsch   if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights when using grid partition\n"));
759566063dSJacob Faibussowitsch   if (targetSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores partition weights when using grid partition\n"));
769566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
779566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
780c569c6eSMark   /* Check grid */
790c569c6eSMark   for (i = 0; i < 3; ++i) Np *= nodes[i] * procs[i];
8063a3b9bcSJacob Faibussowitsch   PetscCheck(nparts == Np, comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %" PetscInt_FMT " grid size", nparts, Np);
8163a3b9bcSJacob Faibussowitsch   PetscCheck(nparts == size, comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %d processes", nparts, size);
8263a3b9bcSJacob Faibussowitsch   PetscCheck(numVertices % nparts == 0, comm, PETSC_ERR_ARG_INCOMP, "Number of cells %" PetscInt_FMT " is not divisible by number of partitions %" PetscInt_FMT, numVertices, nparts);
830c569c6eSMark   for (i = 0; i < p->gridDim; ++i) cells[i] = nodes[i] * procs[i];
840c569c6eSMark   Nr = numVertices / nparts;
850c569c6eSMark   while (Nr > 1) {
860c569c6eSMark     for (i = 0; i < p->gridDim; ++i) {
870c569c6eSMark       cells[i] *= 2;
880c569c6eSMark       Nr /= 2;
890c569c6eSMark     }
900c569c6eSMark   }
911dca8a05SBarry Smith   PetscCheck(!numVertices || Nr == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Odd number of cells %" PetscInt_FMT ". Must be nprocs*2^k", numVertices);
920c569c6eSMark   for (i = 0; i < p->gridDim; ++i) {
9363a3b9bcSJacob Faibussowitsch     PetscCheck(cells[i] % (nodes[i] * procs[i]) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "dir %" PetscInt_FMT ". Number of cells (%" PetscInt_FMT ") mod number of processors %" PetscInt_FMT, i, cells[i], nodes[i] * procs[i]);
94bd026e97SJed Brown     pcells[i] = cells[i] / (nodes[i] * procs[i]);
950c569c6eSMark   }
960c569c6eSMark   /* Compute sizes */
979566063dSJacob Faibussowitsch   for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices / nparts));
989566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(partSection));
999566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(nparts, &offsets));
1009566063dSJacob Faibussowitsch   for (np = 0; np < nparts; ++np) PetscCall(PetscSectionGetOffset(partSection, np, &offsets[np]));
1010c569c6eSMark   if (!numVertices) pcells[0] = pcells[1] = pcells[2] = 0;
1020c569c6eSMark   /* Compute partition */
1039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numVertices, &cellproc));
1040c569c6eSMark   for (nk = 0; nk < nodes[2]; ++nk) {
1050c569c6eSMark     for (nj = 0; nj < nodes[1]; ++nj) {
1060c569c6eSMark       for (ni = 0; ni < nodes[0]; ++ni) {
1070c569c6eSMark         const PetscInt nid = (nk * nodes[1] + nj) * nodes[0] + ni;
1080c569c6eSMark 
1090c569c6eSMark         for (pk = 0; pk < procs[2]; ++pk) {
1100c569c6eSMark           for (pj = 0; pj < procs[1]; ++pj) {
1110c569c6eSMark             for (pi = 0; pi < procs[0]; ++pi) {
1120c569c6eSMark               const PetscInt pid = ((nid * procs[2] + pk) * procs[1] + pj) * procs[0] + pi;
1130c569c6eSMark 
1140c569c6eSMark               /* Assume that cells are originally numbered lexicographically */
1150c569c6eSMark               for (ck = 0; ck < pcells[2]; ++ck) {
1160c569c6eSMark                 for (cj = 0; cj < pcells[1]; ++cj) {
1170c569c6eSMark                   for (ci = 0; ci < pcells[0]; ++ci) {
1180c569c6eSMark                     const PetscInt cid = (((nk * procs[2] + pk) * pcells[2] + ck) * cells[1] + ((nj * procs[1] + pj) * pcells[1] + cj)) * cells[0] + (ni * procs[0] + pi) * pcells[0] + ci;
1190c569c6eSMark 
1200c569c6eSMark                     cellproc[offsets[pid]++] = cid;
1210c569c6eSMark                   }
1220c569c6eSMark                 }
1230c569c6eSMark               }
1240c569c6eSMark             }
1250c569c6eSMark           }
1260c569c6eSMark         }
1270c569c6eSMark       }
1280c569c6eSMark     }
1290c569c6eSMark   }
1301dca8a05SBarry Smith   for (np = 1; np < nparts; ++np) PetscCheck(offsets[np] - offsets[np - 1] == numVertices / nparts, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Offset %" PetscInt_FMT " != %" PetscInt_FMT " partition size", offsets[np], numVertices / nparts);
1319566063dSJacob Faibussowitsch   PetscCall(PetscFree(offsets));
1329566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numVertices, cellproc, PETSC_OWN_POINTER, partition));
133*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1340c569c6eSMark }
1350c569c6eSMark 
136d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
137d71ae5a4SJacob Faibussowitsch {
1380c569c6eSMark   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *)part->data;
139abe9303eSLisandro Dalcin   MPI_Comm                 comm;
140abe9303eSLisandro Dalcin   PetscInt                 np, *tpwgts = NULL, sumw = 0, numVerticesGlobal = 0;
141abe9303eSLisandro Dalcin   PetscMPIInt              size;
142abe9303eSLisandro Dalcin 
143abe9303eSLisandro Dalcin   PetscFunctionBegin;
1440c569c6eSMark   if (p->useGrid) {
1459566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerPartition_Simple_Grid(part, nparts, numVertices, start, adjacency, vertSection, targetSection, partSection, partition));
146*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1470c569c6eSMark   }
1489566063dSJacob Faibussowitsch   if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights\n"));
1499566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
1509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
151abe9303eSLisandro Dalcin   if (targetSection) {
1521c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm));
1539566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(nparts, &tpwgts));
154abe9303eSLisandro Dalcin     for (np = 0; np < nparts; ++np) {
1559566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(targetSection, np, &tpwgts[np]));
156abe9303eSLisandro Dalcin       sumw += tpwgts[np];
157abe9303eSLisandro Dalcin     }
158d7cc930eSLisandro Dalcin     if (sumw) {
159abe9303eSLisandro Dalcin       PetscInt m, mp;
160abe9303eSLisandro Dalcin       for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np] * numVerticesGlobal) / sumw;
161abe9303eSLisandro Dalcin       for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) {
1629371c9d4SSatish Balay         if (m < tpwgts[np]) {
1639371c9d4SSatish Balay           m  = tpwgts[np];
1649371c9d4SSatish Balay           mp = np;
1659371c9d4SSatish Balay         }
166abe9303eSLisandro Dalcin         sumw += tpwgts[np];
167abe9303eSLisandro Dalcin       }
168abe9303eSLisandro Dalcin       if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw;
169abe9303eSLisandro Dalcin     }
1709566063dSJacob Faibussowitsch     if (!sumw) PetscCall(PetscFree(tpwgts));
171abe9303eSLisandro Dalcin   }
172abe9303eSLisandro Dalcin 
1739566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition));
174abe9303eSLisandro Dalcin   if (size == 1) {
175abe9303eSLisandro Dalcin     if (tpwgts) {
17648a46eb9SPierre Jolivet       for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, tpwgts[np]));
177abe9303eSLisandro Dalcin     } else {
17848a46eb9SPierre Jolivet       for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices / nparts + ((numVertices % nparts) > np)));
179abe9303eSLisandro Dalcin     }
180abe9303eSLisandro Dalcin   } else {
181abe9303eSLisandro Dalcin     if (tpwgts) {
182abe9303eSLisandro Dalcin       Vec          v;
183abe9303eSLisandro Dalcin       PetscScalar *array;
184abe9303eSLisandro Dalcin       PetscInt     st, j;
185abe9303eSLisandro Dalcin       PetscMPIInt  rank;
186abe9303eSLisandro Dalcin 
1879566063dSJacob Faibussowitsch       PetscCall(VecCreate(comm, &v));
1889566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(v, numVertices, numVerticesGlobal));
1899566063dSJacob Faibussowitsch       PetscCall(VecSetType(v, VECSTANDARD));
1909566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(comm, &rank));
191abe9303eSLisandro Dalcin       for (np = 0, st = 0; np < nparts; ++np) {
192abe9303eSLisandro Dalcin         if (rank == np || (rank == size - 1 && size < nparts && np >= size)) {
19348a46eb9SPierre Jolivet           for (j = 0; j < tpwgts[np]; j++) PetscCall(VecSetValue(v, st + j, np, INSERT_VALUES));
194abe9303eSLisandro Dalcin         }
195abe9303eSLisandro Dalcin         st += tpwgts[np];
196abe9303eSLisandro Dalcin       }
1979566063dSJacob Faibussowitsch       PetscCall(VecAssemblyBegin(v));
1989566063dSJacob Faibussowitsch       PetscCall(VecAssemblyEnd(v));
1999566063dSJacob Faibussowitsch       PetscCall(VecGetArray(v, &array));
20048a46eb9SPierre Jolivet       for (j = 0; j < numVertices; ++j) PetscCall(PetscSectionAddDof(partSection, PetscRealPart(array[j]), 1));
2019566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(v, &array));
2029566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&v));
203abe9303eSLisandro Dalcin     } else {
204abe9303eSLisandro Dalcin       PetscMPIInt rank;
205abe9303eSLisandro Dalcin       PetscInt    nvGlobal, *offsets, myFirst, myLast;
206abe9303eSLisandro Dalcin 
2079566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(size + 1, &offsets));
208abe9303eSLisandro Dalcin       offsets[0] = 0;
2099566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allgather(&numVertices, 1, MPIU_INT, &offsets[1], 1, MPIU_INT, comm));
210ad540459SPierre Jolivet       for (np = 2; np <= size; np++) offsets[np] += offsets[np - 1];
211abe9303eSLisandro Dalcin       nvGlobal = offsets[size];
2129566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(comm, &rank));
213abe9303eSLisandro Dalcin       myFirst = offsets[rank];
214abe9303eSLisandro Dalcin       myLast  = offsets[rank + 1] - 1;
2159566063dSJacob Faibussowitsch       PetscCall(PetscFree(offsets));
216abe9303eSLisandro Dalcin       if (numVertices) {
217abe9303eSLisandro Dalcin         PetscInt firstPart = 0, firstLargePart = 0;
218abe9303eSLisandro Dalcin         PetscInt lastPart = 0, lastLargePart = 0;
219abe9303eSLisandro Dalcin         PetscInt rem    = nvGlobal % nparts;
220abe9303eSLisandro Dalcin         PetscInt pSmall = nvGlobal / nparts;
221abe9303eSLisandro Dalcin         PetscInt pBig   = nvGlobal / nparts + 1;
222abe9303eSLisandro Dalcin 
223abe9303eSLisandro Dalcin         if (rem) {
224abe9303eSLisandro Dalcin           firstLargePart = myFirst / pBig;
225abe9303eSLisandro Dalcin           lastLargePart  = myLast / pBig;
226abe9303eSLisandro Dalcin 
227abe9303eSLisandro Dalcin           if (firstLargePart < rem) {
228abe9303eSLisandro Dalcin             firstPart = firstLargePart;
229abe9303eSLisandro Dalcin           } else {
230abe9303eSLisandro Dalcin             firstPart = rem + (myFirst - (rem * pBig)) / pSmall;
231abe9303eSLisandro Dalcin           }
232abe9303eSLisandro Dalcin           if (lastLargePart < rem) {
233abe9303eSLisandro Dalcin             lastPart = lastLargePart;
234abe9303eSLisandro Dalcin           } else {
235abe9303eSLisandro Dalcin             lastPart = rem + (myLast - (rem * pBig)) / pSmall;
236abe9303eSLisandro Dalcin           }
237abe9303eSLisandro Dalcin         } else {
238abe9303eSLisandro Dalcin           firstPart = myFirst / (nvGlobal / nparts);
239abe9303eSLisandro Dalcin           lastPart  = myLast / (nvGlobal / nparts);
240abe9303eSLisandro Dalcin         }
241abe9303eSLisandro Dalcin 
242abe9303eSLisandro Dalcin         for (np = firstPart; np <= lastPart; np++) {
243abe9303eSLisandro Dalcin           PetscInt PartStart = np * (nvGlobal / nparts) + PetscMin(nvGlobal % nparts, np);
244abe9303eSLisandro Dalcin           PetscInt PartEnd   = (np + 1) * (nvGlobal / nparts) + PetscMin(nvGlobal % nparts, np + 1);
245abe9303eSLisandro Dalcin 
246abe9303eSLisandro Dalcin           PartStart = PetscMax(PartStart, myFirst);
247abe9303eSLisandro Dalcin           PartEnd   = PetscMin(PartEnd, myLast + 1);
2489566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetDof(partSection, np, PartEnd - PartStart));
249abe9303eSLisandro Dalcin         }
250abe9303eSLisandro Dalcin       }
251abe9303eSLisandro Dalcin     }
252abe9303eSLisandro Dalcin   }
2539566063dSJacob Faibussowitsch   PetscCall(PetscFree(tpwgts));
254*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255abe9303eSLisandro Dalcin }
256abe9303eSLisandro Dalcin 
257d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part)
258d71ae5a4SJacob Faibussowitsch {
259abe9303eSLisandro Dalcin   PetscFunctionBegin;
260abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE;
261abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Simple;
2620c569c6eSMark   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Simple;
263abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Simple;
264abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Simple;
265*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266abe9303eSLisandro Dalcin }
267abe9303eSLisandro Dalcin 
268abe9303eSLisandro Dalcin /*MC
269abe9303eSLisandro Dalcin   PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object
270abe9303eSLisandro Dalcin 
271abe9303eSLisandro Dalcin   Level: intermediate
272abe9303eSLisandro Dalcin 
273db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
274abe9303eSLisandro Dalcin M*/
275abe9303eSLisandro Dalcin 
276d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part)
277d71ae5a4SJacob Faibussowitsch {
278abe9303eSLisandro Dalcin   PetscPartitioner_Simple *p;
279abe9303eSLisandro Dalcin 
280abe9303eSLisandro Dalcin   PetscFunctionBegin;
281abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
2824dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
2830c569c6eSMark   p->gridDim = -1;
284abe9303eSLisandro Dalcin   part->data = p;
285abe9303eSLisandro Dalcin 
2869566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Simple(part));
287*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288abe9303eSLisandro Dalcin }
289