1abe9303eSLisandro Dalcin #include <petscvec.h>
2abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
3abe9303eSLisandro Dalcin
4abe9303eSLisandro Dalcin typedef struct {
50c569c6eSMark PetscBool useGrid; /* Flag to use a grid layout */
60c569c6eSMark PetscInt gridDim; /* The grid dimension */
70c569c6eSMark PetscInt nodeGrid[3]; /* Dimension of node grid */
80c569c6eSMark PetscInt processGrid[3]; /* Dimension of local process grid on each node */
9abe9303eSLisandro Dalcin } PetscPartitioner_Simple;
10abe9303eSLisandro Dalcin
PetscPartitionerDestroy_Simple(PetscPartitioner part)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part)
12d71ae5a4SJacob Faibussowitsch {
13abe9303eSLisandro Dalcin PetscFunctionBegin;
149566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data));
153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16abe9303eSLisandro Dalcin }
17abe9303eSLisandro Dalcin
PetscPartitionerSetFromOptions_Simple(PetscPartitioner part,PetscOptionItems PetscOptionsObject)18*ce78bad3SBarry Smith static PetscErrorCode PetscPartitionerSetFromOptions_Simple(PetscPartitioner part, PetscOptionItems PetscOptionsObject)
19d71ae5a4SJacob Faibussowitsch {
200c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *)part->data;
210c569c6eSMark PetscInt num, i;
220c569c6eSMark PetscBool flg;
230c569c6eSMark
240c569c6eSMark PetscFunctionBegin;
250c569c6eSMark for (i = 0; i < 3; ++i) p->processGrid[i] = p->nodeGrid[i] = 1;
26d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner Simple Options");
270c569c6eSMark num = 3;
289566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_node_grid", "Number of nodes in each dimension", "", p->nodeGrid, &num, &flg));
299371c9d4SSatish Balay if (flg) {
309371c9d4SSatish Balay p->useGrid = PETSC_TRUE;
319371c9d4SSatish Balay p->gridDim = num;
329371c9d4SSatish Balay }
330c569c6eSMark num = 3;
349566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_process_grid", "Number of local processes in each dimension for a given node", "", p->processGrid, &num, &flg));
350c569c6eSMark if (flg) {
360c569c6eSMark p->useGrid = PETSC_TRUE;
370c569c6eSMark if (p->gridDim < 0) p->gridDim = num;
3863a3b9bcSJacob 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);
390c569c6eSMark }
40d0609cedSBarry Smith PetscOptionsHeadEnd();
413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
420c569c6eSMark }
430c569c6eSMark
PetscPartitionerPartition_Simple_Grid(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection targetSection,PetscSection partSection,IS * partition)44d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerPartition_Simple_Grid(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
45d71ae5a4SJacob Faibussowitsch {
460c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *)part->data;
470c569c6eSMark const PetscInt *nodes = p->nodeGrid;
480c569c6eSMark const PetscInt *procs = p->processGrid;
490c569c6eSMark PetscInt *cellproc, *offsets, cells[3] = {1, 1, 1}, pcells[3] = {1, 1, 1};
50bd026e97SJed Brown PetscInt Np = 1, Nr, np, nk, nj, ni, pk, pj, pi, ck, cj, ci, i;
510c569c6eSMark MPI_Comm comm;
520c569c6eSMark PetscMPIInt size;
530c569c6eSMark
540c569c6eSMark PetscFunctionBegin;
559566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights when using grid partition\n"));
569566063dSJacob Faibussowitsch if (targetSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores partition weights when using grid partition\n"));
579566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
589566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
590c569c6eSMark /* Check grid */
600c569c6eSMark for (i = 0; i < 3; ++i) Np *= nodes[i] * procs[i];
6163a3b9bcSJacob Faibussowitsch PetscCheck(nparts == Np, comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %" PetscInt_FMT " grid size", nparts, Np);
6263a3b9bcSJacob Faibussowitsch PetscCheck(nparts == size, comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %d processes", nparts, size);
6363a3b9bcSJacob 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);
640c569c6eSMark for (i = 0; i < p->gridDim; ++i) cells[i] = nodes[i] * procs[i];
650c569c6eSMark Nr = numVertices / nparts;
660c569c6eSMark while (Nr > 1) {
670c569c6eSMark for (i = 0; i < p->gridDim; ++i) {
680c569c6eSMark cells[i] *= 2;
690c569c6eSMark Nr /= 2;
700c569c6eSMark }
710c569c6eSMark }
721dca8a05SBarry Smith PetscCheck(!numVertices || Nr == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Odd number of cells %" PetscInt_FMT ". Must be nprocs*2^k", numVertices);
730c569c6eSMark for (i = 0; i < p->gridDim; ++i) {
7463a3b9bcSJacob 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]);
75bd026e97SJed Brown pcells[i] = cells[i] / (nodes[i] * procs[i]);
760c569c6eSMark }
770c569c6eSMark /* Compute sizes */
789566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices / nparts));
799566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(partSection));
809566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &offsets));
819566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionGetOffset(partSection, np, &offsets[np]));
820c569c6eSMark if (!numVertices) pcells[0] = pcells[1] = pcells[2] = 0;
830c569c6eSMark /* Compute partition */
849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVertices, &cellproc));
850c569c6eSMark for (nk = 0; nk < nodes[2]; ++nk) {
860c569c6eSMark for (nj = 0; nj < nodes[1]; ++nj) {
870c569c6eSMark for (ni = 0; ni < nodes[0]; ++ni) {
880c569c6eSMark const PetscInt nid = (nk * nodes[1] + nj) * nodes[0] + ni;
890c569c6eSMark
900c569c6eSMark for (pk = 0; pk < procs[2]; ++pk) {
910c569c6eSMark for (pj = 0; pj < procs[1]; ++pj) {
920c569c6eSMark for (pi = 0; pi < procs[0]; ++pi) {
930c569c6eSMark const PetscInt pid = ((nid * procs[2] + pk) * procs[1] + pj) * procs[0] + pi;
940c569c6eSMark
950c569c6eSMark /* Assume that cells are originally numbered lexicographically */
960c569c6eSMark for (ck = 0; ck < pcells[2]; ++ck) {
970c569c6eSMark for (cj = 0; cj < pcells[1]; ++cj) {
980c569c6eSMark for (ci = 0; ci < pcells[0]; ++ci) {
990c569c6eSMark 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;
1000c569c6eSMark
1010c569c6eSMark cellproc[offsets[pid]++] = cid;
1020c569c6eSMark }
1030c569c6eSMark }
1040c569c6eSMark }
1050c569c6eSMark }
1060c569c6eSMark }
1070c569c6eSMark }
1080c569c6eSMark }
1090c569c6eSMark }
1100c569c6eSMark }
1111dca8a05SBarry 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);
1129566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets));
1139566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numVertices, cellproc, PETSC_OWN_POINTER, partition));
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1150c569c6eSMark }
1160c569c6eSMark
PetscPartitionerPartition_Simple(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * partition)11721c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
118d71ae5a4SJacob Faibussowitsch {
1190c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *)part->data;
120abe9303eSLisandro Dalcin MPI_Comm comm;
121abe9303eSLisandro Dalcin PetscInt np, *tpwgts = NULL, sumw = 0, numVerticesGlobal = 0;
122abe9303eSLisandro Dalcin PetscMPIInt size;
123abe9303eSLisandro Dalcin
124abe9303eSLisandro Dalcin PetscFunctionBegin;
1250c569c6eSMark if (p->useGrid) {
1269566063dSJacob Faibussowitsch PetscCall(PetscPartitionerPartition_Simple_Grid(part, nparts, numVertices, start, adjacency, vertSection, targetSection, partSection, partition));
1273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
1280c569c6eSMark }
1299566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights\n"));
1309566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
132abe9303eSLisandro Dalcin if (targetSection) {
133462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm));
1349566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &tpwgts));
135abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) {
1369566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection, np, &tpwgts[np]));
137abe9303eSLisandro Dalcin sumw += tpwgts[np];
138abe9303eSLisandro Dalcin }
139d7cc930eSLisandro Dalcin if (sumw) {
140abe9303eSLisandro Dalcin PetscInt m, mp;
141abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np] * numVerticesGlobal) / sumw;
142abe9303eSLisandro Dalcin for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) {
1439371c9d4SSatish Balay if (m < tpwgts[np]) {
1449371c9d4SSatish Balay m = tpwgts[np];
1459371c9d4SSatish Balay mp = np;
1469371c9d4SSatish Balay }
147abe9303eSLisandro Dalcin sumw += tpwgts[np];
148abe9303eSLisandro Dalcin }
149abe9303eSLisandro Dalcin if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw;
150abe9303eSLisandro Dalcin }
1519566063dSJacob Faibussowitsch if (!sumw) PetscCall(PetscFree(tpwgts));
152abe9303eSLisandro Dalcin }
153abe9303eSLisandro Dalcin
1549566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition));
155abe9303eSLisandro Dalcin if (size == 1) {
156abe9303eSLisandro Dalcin if (tpwgts) {
15748a46eb9SPierre Jolivet for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, tpwgts[np]));
158abe9303eSLisandro Dalcin } else {
15948a46eb9SPierre Jolivet for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices / nparts + ((numVertices % nparts) > np)));
160abe9303eSLisandro Dalcin }
161abe9303eSLisandro Dalcin } else {
162abe9303eSLisandro Dalcin if (tpwgts) {
163abe9303eSLisandro Dalcin Vec v;
164abe9303eSLisandro Dalcin PetscScalar *array;
165abe9303eSLisandro Dalcin PetscInt st, j;
166abe9303eSLisandro Dalcin PetscMPIInt rank;
167abe9303eSLisandro Dalcin
1689566063dSJacob Faibussowitsch PetscCall(VecCreate(comm, &v));
1699566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v, numVertices, numVerticesGlobal));
1709566063dSJacob Faibussowitsch PetscCall(VecSetType(v, VECSTANDARD));
1719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
172abe9303eSLisandro Dalcin for (np = 0, st = 0; np < nparts; ++np) {
173abe9303eSLisandro Dalcin if (rank == np || (rank == size - 1 && size < nparts && np >= size)) {
17448a46eb9SPierre Jolivet for (j = 0; j < tpwgts[np]; j++) PetscCall(VecSetValue(v, st + j, np, INSERT_VALUES));
175abe9303eSLisandro Dalcin }
176abe9303eSLisandro Dalcin st += tpwgts[np];
177abe9303eSLisandro Dalcin }
1789566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(v));
1799566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(v));
1809566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &array));
18148a46eb9SPierre Jolivet for (j = 0; j < numVertices; ++j) PetscCall(PetscSectionAddDof(partSection, PetscRealPart(array[j]), 1));
1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &array));
1839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v));
184abe9303eSLisandro Dalcin } else {
185abe9303eSLisandro Dalcin PetscMPIInt rank;
186abe9303eSLisandro Dalcin PetscInt nvGlobal, *offsets, myFirst, myLast;
187abe9303eSLisandro Dalcin
1889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size + 1, &offsets));
189abe9303eSLisandro Dalcin offsets[0] = 0;
1909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&numVertices, 1, MPIU_INT, &offsets[1], 1, MPIU_INT, comm));
191ad540459SPierre Jolivet for (np = 2; np <= size; np++) offsets[np] += offsets[np - 1];
192abe9303eSLisandro Dalcin nvGlobal = offsets[size];
1939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
194abe9303eSLisandro Dalcin myFirst = offsets[rank];
195abe9303eSLisandro Dalcin myLast = offsets[rank + 1] - 1;
1969566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets));
197abe9303eSLisandro Dalcin if (numVertices) {
198abe9303eSLisandro Dalcin PetscInt firstPart = 0, firstLargePart = 0;
199abe9303eSLisandro Dalcin PetscInt lastPart = 0, lastLargePart = 0;
200abe9303eSLisandro Dalcin PetscInt rem = nvGlobal % nparts;
201abe9303eSLisandro Dalcin PetscInt pSmall = nvGlobal / nparts;
202abe9303eSLisandro Dalcin PetscInt pBig = nvGlobal / nparts + 1;
203abe9303eSLisandro Dalcin
204abe9303eSLisandro Dalcin if (rem) {
205abe9303eSLisandro Dalcin firstLargePart = myFirst / pBig;
206abe9303eSLisandro Dalcin lastLargePart = myLast / pBig;
207abe9303eSLisandro Dalcin
208abe9303eSLisandro Dalcin if (firstLargePart < rem) {
209abe9303eSLisandro Dalcin firstPart = firstLargePart;
210abe9303eSLisandro Dalcin } else {
211abe9303eSLisandro Dalcin firstPart = rem + (myFirst - (rem * pBig)) / pSmall;
212abe9303eSLisandro Dalcin }
213abe9303eSLisandro Dalcin if (lastLargePart < rem) {
214abe9303eSLisandro Dalcin lastPart = lastLargePart;
215abe9303eSLisandro Dalcin } else {
216abe9303eSLisandro Dalcin lastPart = rem + (myLast - (rem * pBig)) / pSmall;
217abe9303eSLisandro Dalcin }
218abe9303eSLisandro Dalcin } else {
219abe9303eSLisandro Dalcin firstPart = myFirst / (nvGlobal / nparts);
220abe9303eSLisandro Dalcin lastPart = myLast / (nvGlobal / nparts);
221abe9303eSLisandro Dalcin }
222abe9303eSLisandro Dalcin
223abe9303eSLisandro Dalcin for (np = firstPart; np <= lastPart; np++) {
224abe9303eSLisandro Dalcin PetscInt PartStart = np * (nvGlobal / nparts) + PetscMin(nvGlobal % nparts, np);
225abe9303eSLisandro Dalcin PetscInt PartEnd = (np + 1) * (nvGlobal / nparts) + PetscMin(nvGlobal % nparts, np + 1);
226abe9303eSLisandro Dalcin
227abe9303eSLisandro Dalcin PartStart = PetscMax(PartStart, myFirst);
228abe9303eSLisandro Dalcin PartEnd = PetscMin(PartEnd, myLast + 1);
2299566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, np, PartEnd - PartStart));
230abe9303eSLisandro Dalcin }
231abe9303eSLisandro Dalcin }
232abe9303eSLisandro Dalcin }
233abe9303eSLisandro Dalcin }
2349566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts));
2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
236abe9303eSLisandro Dalcin }
237abe9303eSLisandro Dalcin
PetscPartitionerInitialize_Simple(PetscPartitioner part)238d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part)
239d71ae5a4SJacob Faibussowitsch {
240abe9303eSLisandro Dalcin PetscFunctionBegin;
241abe9303eSLisandro Dalcin part->noGraph = PETSC_TRUE;
2420c569c6eSMark part->ops->setfromoptions = PetscPartitionerSetFromOptions_Simple;
243abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_Simple;
244abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_Simple;
2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
246abe9303eSLisandro Dalcin }
247abe9303eSLisandro Dalcin
248abe9303eSLisandro Dalcin /*MC
249abe9303eSLisandro Dalcin PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object
250abe9303eSLisandro Dalcin
251abe9303eSLisandro Dalcin Level: intermediate
252abe9303eSLisandro Dalcin
253db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
254abe9303eSLisandro Dalcin M*/
255abe9303eSLisandro Dalcin
PetscPartitionerCreate_Simple(PetscPartitioner part)256d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part)
257d71ae5a4SJacob Faibussowitsch {
258abe9303eSLisandro Dalcin PetscPartitioner_Simple *p;
259abe9303eSLisandro Dalcin
260abe9303eSLisandro Dalcin PetscFunctionBegin;
261abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
2624dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&p));
2630c569c6eSMark p->gridDim = -1;
264abe9303eSLisandro Dalcin part->data = p;
265abe9303eSLisandro Dalcin
2669566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_Simple(part));
2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
268abe9303eSLisandro Dalcin }
269