xref: /petsc/src/dm/partitioner/impls/simple/partsimple.c (revision 0c569c6e26eeef44c3acfb80fe8a62dba0a07e6c)
1abe9303eSLisandro Dalcin #include <petscvec.h>
2abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h>        /*I "petscpartitioner.h" I*/
3abe9303eSLisandro Dalcin 
4abe9303eSLisandro Dalcin typedef struct {
5*0c569c6eSMark   PetscBool useGrid;        /* Flag to use a grid layout */
6*0c569c6eSMark   PetscInt  gridDim;        /* The grid dimension */
7*0c569c6eSMark   PetscInt  nodeGrid[3];    /* Dimension of node grid */
8*0c569c6eSMark   PetscInt  processGrid[3]; /* Dimension of local process grid on each node */
9abe9303eSLisandro Dalcin } PetscPartitioner_Simple;
10abe9303eSLisandro Dalcin 
11abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part)
12abe9303eSLisandro Dalcin {
13abe9303eSLisandro Dalcin   PetscErrorCode ierr;
14abe9303eSLisandro Dalcin 
15abe9303eSLisandro Dalcin   PetscFunctionBegin;
16abe9303eSLisandro Dalcin   ierr = PetscFree(part->data);CHKERRQ(ierr);
17abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
18abe9303eSLisandro Dalcin }
19abe9303eSLisandro Dalcin 
20abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple_ASCII(PetscPartitioner part, PetscViewer viewer)
21abe9303eSLisandro Dalcin {
22abe9303eSLisandro Dalcin   PetscFunctionBegin;
23abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
24abe9303eSLisandro Dalcin }
25abe9303eSLisandro Dalcin 
26abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer)
27abe9303eSLisandro Dalcin {
28abe9303eSLisandro Dalcin   PetscBool      iascii;
29abe9303eSLisandro Dalcin   PetscErrorCode ierr;
30abe9303eSLisandro Dalcin 
31abe9303eSLisandro Dalcin   PetscFunctionBegin;
32abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
33abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
34abe9303eSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
35abe9303eSLisandro Dalcin   if (iascii) {ierr = PetscPartitionerView_Simple_ASCII(part, viewer);CHKERRQ(ierr);}
36abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
37abe9303eSLisandro Dalcin }
38abe9303eSLisandro Dalcin 
39*0c569c6eSMark static PetscErrorCode PetscPartitionerSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
40*0c569c6eSMark {
41*0c569c6eSMark   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data;
42*0c569c6eSMark   PetscInt                 num, i;
43*0c569c6eSMark   PetscBool                flg;
44*0c569c6eSMark   PetscErrorCode           ierr;
45*0c569c6eSMark 
46*0c569c6eSMark   PetscFunctionBegin;
47*0c569c6eSMark   for (i = 0; i < 3; ++i) p->processGrid[i] = p->nodeGrid[i] = 1;
48*0c569c6eSMark   ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitioner Simple Options");CHKERRQ(ierr);
49*0c569c6eSMark   num  = 3;
50*0c569c6eSMark   ierr = PetscOptionsIntArray("-petscpartitioner_simple_node_grid", "Number of nodes in each dimension", "", p->nodeGrid, &num, &flg);CHKERRQ(ierr);
51*0c569c6eSMark   if (flg) {p->useGrid = PETSC_TRUE; p->gridDim = num;}
52*0c569c6eSMark   num  = 3;
53*0c569c6eSMark   ierr = PetscOptionsIntArray("-petscpartitioner_simple_process_grid", "Number of local processes in each dimension for a given node", "", p->processGrid, &num, &flg);CHKERRQ(ierr);
54*0c569c6eSMark   if (flg) {
55*0c569c6eSMark     p->useGrid = PETSC_TRUE;
56*0c569c6eSMark     if (p->gridDim < 0) p->gridDim = num;
57*0c569c6eSMark     else if (p->gridDim != num) SETERRQ2(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_INCOMP, "Process grid dimension %D != %D node grid dimension", num, p->gridDim);
58*0c569c6eSMark   }
59*0c569c6eSMark   ierr = PetscOptionsTail();CHKERRQ(ierr);
60*0c569c6eSMark   PetscFunctionReturn(0);
61*0c569c6eSMark }
62*0c569c6eSMark 
63*0c569c6eSMark static PetscErrorCode PetscPartitionerPartition_Simple_Grid(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
64*0c569c6eSMark {
65*0c569c6eSMark   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data;
66*0c569c6eSMark   const PetscInt          *nodes = p->nodeGrid;
67*0c569c6eSMark   const PetscInt          *procs = p->processGrid;
68*0c569c6eSMark   PetscInt                *cellproc, *offsets, cells[3] = {1, 1, 1}, pcells[3] = {1, 1, 1};
69*0c569c6eSMark   PetscInt                 Np    = 1, Nlc = 1, Nr, np, nk, nj, ni, pk, pj, pi, ck, cj, ci, i;
70*0c569c6eSMark   MPI_Comm                 comm;
71*0c569c6eSMark   PetscMPIInt              size;
72*0c569c6eSMark   PetscErrorCode           ierr;
73*0c569c6eSMark 
74*0c569c6eSMark   PetscFunctionBegin;
75*0c569c6eSMark   if (vertSection)   {ierr = PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights when using grid partition\n");CHKERRQ(ierr);}
76*0c569c6eSMark   if (targetSection) {ierr = PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores partition weights when using grid partition\n");CHKERRQ(ierr);}
77*0c569c6eSMark   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
78*0c569c6eSMark   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
79*0c569c6eSMark   /* Check grid */
80*0c569c6eSMark   for (i = 0; i < 3; ++i) Np *= nodes[i]*procs[i];
81*0c569c6eSMark   if (nparts != Np)   SETERRQ2(comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %D != %D grid size", nparts, Np);
82*0c569c6eSMark   if (nparts != size) SETERRQ2(comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %D != %D processes", nparts, size);
83*0c569c6eSMark   if (numVertices % nparts) SETERRQ2(comm, PETSC_ERR_ARG_INCOMP, "Number of cells %D is not divisible by number of partitions %D", numVertices, nparts);
84*0c569c6eSMark   for (i = 0; i < p->gridDim; ++i) cells[i] = nodes[i]*procs[i];
85*0c569c6eSMark   Nr = numVertices / nparts;
86*0c569c6eSMark   while (Nr > 1) {
87*0c569c6eSMark     for (i = 0; i < p->gridDim; ++i) {
88*0c569c6eSMark       cells[i] *= 2;
89*0c569c6eSMark       Nr       /= 2;
90*0c569c6eSMark     }
91*0c569c6eSMark   }
92*0c569c6eSMark   if (numVertices && Nr != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Odd number of cells %D. Must be nprocs*2^k", numVertices);
93*0c569c6eSMark   for (i = 0; i < p->gridDim; ++i) {
94*0c569c6eSMark     if (cells[i] %  (nodes[i]*procs[i])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "dir %D. Number of cells (%D) mod number of processors %D", i, cells[i], nodes[i]*procs[i]);
95*0c569c6eSMark     pcells[i] = cells[i] / (nodes[i]*procs[i]); Nlc *= pcells[i];
96*0c569c6eSMark   }
97*0c569c6eSMark   /* Compute sizes */
98*0c569c6eSMark   for (np = 0; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, numVertices/nparts);CHKERRQ(ierr);}
99*0c569c6eSMark   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
100*0c569c6eSMark   ierr = PetscCalloc1(nparts, &offsets);CHKERRQ(ierr);
101*0c569c6eSMark   for (np = 0; np < nparts; ++np) {ierr = PetscSectionGetOffset(partSection, np, &offsets[np]);CHKERRQ(ierr);}
102*0c569c6eSMark   if (!numVertices) pcells[0] = pcells[1] = pcells[2] = 0;
103*0c569c6eSMark   /* Compute partition */
104*0c569c6eSMark   ierr = PetscMalloc1(numVertices, &cellproc);CHKERRQ(ierr);
105*0c569c6eSMark   for (nk = 0; nk < nodes[2]; ++nk) {
106*0c569c6eSMark     for (nj = 0; nj < nodes[1]; ++nj) {
107*0c569c6eSMark       for (ni = 0; ni < nodes[0]; ++ni) {
108*0c569c6eSMark         const PetscInt nid = (nk*nodes[1] + nj)*nodes[0] + ni;
109*0c569c6eSMark 
110*0c569c6eSMark         for (pk = 0; pk < procs[2]; ++pk) {
111*0c569c6eSMark           for (pj = 0; pj < procs[1]; ++pj) {
112*0c569c6eSMark             for (pi = 0; pi < procs[0]; ++pi) {
113*0c569c6eSMark               const PetscInt pid = ((nid*procs[2] + pk)*procs[1] + pj)*procs[0] + pi;
114*0c569c6eSMark 
115*0c569c6eSMark               /* Assume that cells are originally numbered lexicographically */
116*0c569c6eSMark               for (ck = 0; ck < pcells[2]; ++ck) {
117*0c569c6eSMark                 for (cj = 0; cj < pcells[1]; ++cj) {
118*0c569c6eSMark                   for (ci = 0; ci < pcells[0]; ++ci) {
119*0c569c6eSMark                     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;
120*0c569c6eSMark 
121*0c569c6eSMark                     cellproc[offsets[pid]++] = cid;
122*0c569c6eSMark                   }
123*0c569c6eSMark                 }
124*0c569c6eSMark               }
125*0c569c6eSMark             }
126*0c569c6eSMark           }
127*0c569c6eSMark         }
128*0c569c6eSMark       }
129*0c569c6eSMark     }
130*0c569c6eSMark   }
131*0c569c6eSMark   for (np = 1; np < nparts; ++np) if (offsets[np] - offsets[np-1] != numVertices/nparts) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Offset %D != %D partition size", offsets[np], numVertices/nparts);
132*0c569c6eSMark   ierr = PetscFree(offsets);CHKERRQ(ierr);
133*0c569c6eSMark   ierr = ISCreateGeneral(PETSC_COMM_SELF, numVertices, cellproc, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
134*0c569c6eSMark   PetscFunctionReturn(0);
135*0c569c6eSMark }
136*0c569c6eSMark 
137abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
138abe9303eSLisandro Dalcin {
139*0c569c6eSMark   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data;
140abe9303eSLisandro Dalcin   MPI_Comm       comm;
141abe9303eSLisandro Dalcin   PetscInt       np, *tpwgts = NULL, sumw = 0, numVerticesGlobal  = 0;
142abe9303eSLisandro Dalcin   PetscMPIInt    size;
143abe9303eSLisandro Dalcin   PetscErrorCode ierr;
144abe9303eSLisandro Dalcin 
145abe9303eSLisandro Dalcin   PetscFunctionBegin;
146*0c569c6eSMark   if (p->useGrid) {
147*0c569c6eSMark     ierr = PetscPartitionerPartition_Simple_Grid(part, nparts, numVertices, start, adjacency, vertSection, targetSection, partSection, partition);CHKERRQ(ierr);
148*0c569c6eSMark     PetscFunctionReturn(0);
149*0c569c6eSMark   }
150abe9303eSLisandro Dalcin   if (vertSection) {ierr = PetscInfo(part,"PETSCPARTITIONERSIMPLE ignores vertex weights\n");CHKERRQ(ierr);}
151*0c569c6eSMark   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
152abe9303eSLisandro Dalcin   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
153abe9303eSLisandro Dalcin   if (targetSection) {
154abe9303eSLisandro Dalcin     ierr = MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr);
155abe9303eSLisandro Dalcin     ierr = PetscCalloc1(nparts,&tpwgts);CHKERRQ(ierr);
156abe9303eSLisandro Dalcin     for (np = 0; np < nparts; ++np) {
157abe9303eSLisandro Dalcin       ierr = PetscSectionGetDof(targetSection,np,&tpwgts[np]);CHKERRQ(ierr);
158abe9303eSLisandro Dalcin       sumw += tpwgts[np];
159abe9303eSLisandro Dalcin     }
160d7cc930eSLisandro Dalcin     if (sumw) {
161abe9303eSLisandro Dalcin       PetscInt m,mp;
162abe9303eSLisandro Dalcin       for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np]*numVerticesGlobal)/sumw;
163abe9303eSLisandro Dalcin       for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) {
164abe9303eSLisandro Dalcin         if (m < tpwgts[np]) { m = tpwgts[np]; mp = np; }
165abe9303eSLisandro Dalcin         sumw += tpwgts[np];
166abe9303eSLisandro Dalcin       }
167abe9303eSLisandro Dalcin       if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw;
168abe9303eSLisandro Dalcin     }
169d7cc930eSLisandro Dalcin     if (!sumw) {ierr = PetscFree(tpwgts);CHKERRQ(ierr);}
170abe9303eSLisandro Dalcin   }
171abe9303eSLisandro Dalcin 
172abe9303eSLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
173abe9303eSLisandro Dalcin   if (size == 1) {
174abe9303eSLisandro Dalcin     if (tpwgts) {
175abe9303eSLisandro Dalcin       for (np = 0; np < nparts; ++np) {
176abe9303eSLisandro Dalcin         ierr = PetscSectionSetDof(partSection, np, tpwgts[np]);CHKERRQ(ierr);
177abe9303eSLisandro Dalcin       }
178abe9303eSLisandro Dalcin     } else {
179abe9303eSLisandro Dalcin       for (np = 0; np < nparts; ++np) {
180abe9303eSLisandro Dalcin         ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr);
181abe9303eSLisandro Dalcin       }
182abe9303eSLisandro Dalcin     }
183abe9303eSLisandro Dalcin   } else {
184abe9303eSLisandro Dalcin     if (tpwgts) {
185abe9303eSLisandro Dalcin       Vec         v;
186abe9303eSLisandro Dalcin       PetscScalar *array;
187abe9303eSLisandro Dalcin       PetscInt    st,j;
188abe9303eSLisandro Dalcin       PetscMPIInt rank;
189abe9303eSLisandro Dalcin 
190abe9303eSLisandro Dalcin       ierr = VecCreate(comm,&v);CHKERRQ(ierr);
191abe9303eSLisandro Dalcin       ierr = VecSetSizes(v,numVertices,numVerticesGlobal);CHKERRQ(ierr);
192abe9303eSLisandro Dalcin       ierr = VecSetType(v,VECSTANDARD);CHKERRQ(ierr);
193abe9303eSLisandro Dalcin       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
194abe9303eSLisandro Dalcin       for (np = 0,st = 0; np < nparts; ++np) {
195abe9303eSLisandro Dalcin         if (rank == np || (rank == size-1 && size < nparts && np >= size)) {
196abe9303eSLisandro Dalcin           for (j = 0; j < tpwgts[np]; j++) {
197abe9303eSLisandro Dalcin             ierr = VecSetValue(v,st+j,np,INSERT_VALUES);CHKERRQ(ierr);
198abe9303eSLisandro Dalcin           }
199abe9303eSLisandro Dalcin         }
200abe9303eSLisandro Dalcin         st += tpwgts[np];
201abe9303eSLisandro Dalcin       }
202abe9303eSLisandro Dalcin       ierr = VecAssemblyBegin(v);CHKERRQ(ierr);
203abe9303eSLisandro Dalcin       ierr = VecAssemblyEnd(v);CHKERRQ(ierr);
204abe9303eSLisandro Dalcin       ierr = VecGetArray(v,&array);CHKERRQ(ierr);
205abe9303eSLisandro Dalcin       for (j = 0; j < numVertices; ++j) {
206abe9303eSLisandro Dalcin         ierr = PetscSectionAddDof(partSection,PetscRealPart(array[j]),1);CHKERRQ(ierr);
207abe9303eSLisandro Dalcin       }
208abe9303eSLisandro Dalcin       ierr = VecRestoreArray(v,&array);CHKERRQ(ierr);
209abe9303eSLisandro Dalcin       ierr = VecDestroy(&v);CHKERRQ(ierr);
210abe9303eSLisandro Dalcin     } else {
211abe9303eSLisandro Dalcin       PetscMPIInt rank;
212abe9303eSLisandro Dalcin       PetscInt nvGlobal, *offsets, myFirst, myLast;
213abe9303eSLisandro Dalcin 
214abe9303eSLisandro Dalcin       ierr = PetscMalloc1(size+1,&offsets);CHKERRQ(ierr);
215abe9303eSLisandro Dalcin       offsets[0] = 0;
216abe9303eSLisandro Dalcin       ierr = MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm);CHKERRQ(ierr);
217abe9303eSLisandro Dalcin       for (np = 2; np <= size; np++) {
218abe9303eSLisandro Dalcin         offsets[np] += offsets[np-1];
219abe9303eSLisandro Dalcin       }
220abe9303eSLisandro Dalcin       nvGlobal = offsets[size];
221abe9303eSLisandro Dalcin       ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
222abe9303eSLisandro Dalcin       myFirst = offsets[rank];
223abe9303eSLisandro Dalcin       myLast  = offsets[rank + 1] - 1;
224abe9303eSLisandro Dalcin       ierr = PetscFree(offsets);CHKERRQ(ierr);
225abe9303eSLisandro Dalcin       if (numVertices) {
226abe9303eSLisandro Dalcin         PetscInt firstPart = 0, firstLargePart = 0;
227abe9303eSLisandro Dalcin         PetscInt lastPart = 0, lastLargePart = 0;
228abe9303eSLisandro Dalcin         PetscInt rem = nvGlobal % nparts;
229abe9303eSLisandro Dalcin         PetscInt pSmall = nvGlobal/nparts;
230abe9303eSLisandro Dalcin         PetscInt pBig = nvGlobal/nparts + 1;
231abe9303eSLisandro Dalcin 
232abe9303eSLisandro Dalcin         if (rem) {
233abe9303eSLisandro Dalcin           firstLargePart = myFirst / pBig;
234abe9303eSLisandro Dalcin           lastLargePart  = myLast  / pBig;
235abe9303eSLisandro Dalcin 
236abe9303eSLisandro Dalcin           if (firstLargePart < rem) {
237abe9303eSLisandro Dalcin             firstPart = firstLargePart;
238abe9303eSLisandro Dalcin           } else {
239abe9303eSLisandro Dalcin             firstPart = rem + (myFirst - (rem * pBig)) / pSmall;
240abe9303eSLisandro Dalcin           }
241abe9303eSLisandro Dalcin           if (lastLargePart < rem) {
242abe9303eSLisandro Dalcin             lastPart = lastLargePart;
243abe9303eSLisandro Dalcin           } else {
244abe9303eSLisandro Dalcin             lastPart = rem + (myLast - (rem * pBig)) / pSmall;
245abe9303eSLisandro Dalcin           }
246abe9303eSLisandro Dalcin         } else {
247abe9303eSLisandro Dalcin           firstPart = myFirst / (nvGlobal/nparts);
248abe9303eSLisandro Dalcin           lastPart  = myLast  / (nvGlobal/nparts);
249abe9303eSLisandro Dalcin         }
250abe9303eSLisandro Dalcin 
251abe9303eSLisandro Dalcin         for (np = firstPart; np <= lastPart; np++) {
252abe9303eSLisandro Dalcin           PetscInt PartStart =  np    * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np);
253abe9303eSLisandro Dalcin           PetscInt PartEnd   = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1);
254abe9303eSLisandro Dalcin 
255abe9303eSLisandro Dalcin           PartStart = PetscMax(PartStart,myFirst);
256abe9303eSLisandro Dalcin           PartEnd   = PetscMin(PartEnd,myLast+1);
257abe9303eSLisandro Dalcin           ierr = PetscSectionSetDof(partSection,np,PartEnd-PartStart);CHKERRQ(ierr);
258abe9303eSLisandro Dalcin         }
259abe9303eSLisandro Dalcin       }
260abe9303eSLisandro Dalcin     }
261abe9303eSLisandro Dalcin   }
262abe9303eSLisandro Dalcin   ierr = PetscFree(tpwgts);CHKERRQ(ierr);
263abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
264abe9303eSLisandro Dalcin }
265abe9303eSLisandro Dalcin 
266abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part)
267abe9303eSLisandro Dalcin {
268abe9303eSLisandro Dalcin   PetscFunctionBegin;
269abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE;
270abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Simple;
271*0c569c6eSMark   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Simple;
272abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Simple;
273abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Simple;
274abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
275abe9303eSLisandro Dalcin }
276abe9303eSLisandro Dalcin 
277abe9303eSLisandro Dalcin /*MC
278abe9303eSLisandro Dalcin   PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object
279abe9303eSLisandro Dalcin 
280abe9303eSLisandro Dalcin   Level: intermediate
281abe9303eSLisandro Dalcin 
282abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
283abe9303eSLisandro Dalcin M*/
284abe9303eSLisandro Dalcin 
285abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part)
286abe9303eSLisandro Dalcin {
287abe9303eSLisandro Dalcin   PetscPartitioner_Simple *p;
288abe9303eSLisandro Dalcin   PetscErrorCode           ierr;
289abe9303eSLisandro Dalcin 
290abe9303eSLisandro Dalcin   PetscFunctionBegin;
291abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
292abe9303eSLisandro Dalcin   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
293*0c569c6eSMark   p->gridDim = -1;
294abe9303eSLisandro Dalcin   part->data = p;
295abe9303eSLisandro Dalcin 
296abe9303eSLisandro Dalcin   ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr);
297abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
298abe9303eSLisandro Dalcin }
299