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 11abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part) 12abe9303eSLisandro Dalcin { 13abe9303eSLisandro Dalcin PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 15abe9303eSLisandro Dalcin PetscFunctionReturn(0); 16abe9303eSLisandro Dalcin } 17abe9303eSLisandro Dalcin 18abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple_ASCII(PetscPartitioner part, PetscViewer viewer) 19abe9303eSLisandro Dalcin { 20abe9303eSLisandro Dalcin PetscFunctionBegin; 21abe9303eSLisandro Dalcin PetscFunctionReturn(0); 22abe9303eSLisandro Dalcin } 23abe9303eSLisandro Dalcin 24abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer) 25abe9303eSLisandro Dalcin { 26abe9303eSLisandro Dalcin PetscBool iascii; 27abe9303eSLisandro Dalcin 28abe9303eSLisandro Dalcin PetscFunctionBegin; 29abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 30abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 329566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscPartitionerView_Simple_ASCII(part, viewer)); 33abe9303eSLisandro Dalcin PetscFunctionReturn(0); 34abe9303eSLisandro Dalcin } 35abe9303eSLisandro Dalcin 360c569c6eSMark static PetscErrorCode PetscPartitionerSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject, PetscPartitioner part) 370c569c6eSMark { 380c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 390c569c6eSMark PetscInt num, i; 400c569c6eSMark PetscBool flg; 410c569c6eSMark 420c569c6eSMark PetscFunctionBegin; 430c569c6eSMark for (i = 0; i < 3; ++i) p->processGrid[i] = p->nodeGrid[i] = 1; 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsHead(PetscOptionsObject, "PetscPartitioner Simple Options")); 450c569c6eSMark num = 3; 469566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_node_grid", "Number of nodes in each dimension", "", p->nodeGrid, &num, &flg)); 470c569c6eSMark if (flg) {p->useGrid = PETSC_TRUE; p->gridDim = num;} 480c569c6eSMark num = 3; 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_process_grid", "Number of local processes in each dimension for a given node", "", p->processGrid, &num, &flg)); 500c569c6eSMark if (flg) { 510c569c6eSMark p->useGrid = PETSC_TRUE; 520c569c6eSMark if (p->gridDim < 0) p->gridDim = num; 532c71b3e2SJacob Faibussowitsch else PetscCheckFalse(p->gridDim != num,PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_INCOMP, "Process grid dimension %D != %D node grid dimension", num, p->gridDim); 540c569c6eSMark } 559566063dSJacob Faibussowitsch PetscCall(PetscOptionsTail()); 560c569c6eSMark PetscFunctionReturn(0); 570c569c6eSMark } 580c569c6eSMark 590c569c6eSMark static PetscErrorCode PetscPartitionerPartition_Simple_Grid(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 600c569c6eSMark { 610c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 620c569c6eSMark const PetscInt *nodes = p->nodeGrid; 630c569c6eSMark const PetscInt *procs = p->processGrid; 640c569c6eSMark PetscInt *cellproc, *offsets, cells[3] = {1, 1, 1}, pcells[3] = {1, 1, 1}; 65bd026e97SJed Brown PetscInt Np = 1, Nr, np, nk, nj, ni, pk, pj, pi, ck, cj, ci, i; 660c569c6eSMark MPI_Comm comm; 670c569c6eSMark PetscMPIInt size; 680c569c6eSMark 690c569c6eSMark PetscFunctionBegin; 709566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights when using grid partition\n")); 719566063dSJacob Faibussowitsch if (targetSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores partition weights when using grid partition\n")); 729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) part, &comm)); 739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 740c569c6eSMark /* Check grid */ 750c569c6eSMark for (i = 0; i < 3; ++i) Np *= nodes[i]*procs[i]; 762c71b3e2SJacob Faibussowitsch PetscCheckFalse(nparts != Np,comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %D != %D grid size", nparts, Np); 772c71b3e2SJacob Faibussowitsch PetscCheckFalse(nparts != size,comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %D != %D processes", nparts, size); 782c71b3e2SJacob Faibussowitsch PetscCheckFalse(numVertices % nparts,comm, PETSC_ERR_ARG_INCOMP, "Number of cells %D is not divisible by number of partitions %D", numVertices, nparts); 790c569c6eSMark for (i = 0; i < p->gridDim; ++i) cells[i] = nodes[i]*procs[i]; 800c569c6eSMark Nr = numVertices / nparts; 810c569c6eSMark while (Nr > 1) { 820c569c6eSMark for (i = 0; i < p->gridDim; ++i) { 830c569c6eSMark cells[i] *= 2; 840c569c6eSMark Nr /= 2; 850c569c6eSMark } 860c569c6eSMark } 872c71b3e2SJacob Faibussowitsch PetscCheckFalse(numVertices && Nr != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Odd number of cells %D. Must be nprocs*2^k", numVertices); 880c569c6eSMark for (i = 0; i < p->gridDim; ++i) { 892c71b3e2SJacob Faibussowitsch PetscCheckFalse(cells[i] % (nodes[i]*procs[i]),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]); 90bd026e97SJed Brown pcells[i] = cells[i] / (nodes[i]*procs[i]); 910c569c6eSMark } 920c569c6eSMark /* Compute sizes */ 939566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices/nparts)); 949566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(partSection)); 959566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &offsets)); 969566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionGetOffset(partSection, np, &offsets[np])); 970c569c6eSMark if (!numVertices) pcells[0] = pcells[1] = pcells[2] = 0; 980c569c6eSMark /* Compute partition */ 999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVertices, &cellproc)); 1000c569c6eSMark for (nk = 0; nk < nodes[2]; ++nk) { 1010c569c6eSMark for (nj = 0; nj < nodes[1]; ++nj) { 1020c569c6eSMark for (ni = 0; ni < nodes[0]; ++ni) { 1030c569c6eSMark const PetscInt nid = (nk*nodes[1] + nj)*nodes[0] + ni; 1040c569c6eSMark 1050c569c6eSMark for (pk = 0; pk < procs[2]; ++pk) { 1060c569c6eSMark for (pj = 0; pj < procs[1]; ++pj) { 1070c569c6eSMark for (pi = 0; pi < procs[0]; ++pi) { 1080c569c6eSMark const PetscInt pid = ((nid*procs[2] + pk)*procs[1] + pj)*procs[0] + pi; 1090c569c6eSMark 1100c569c6eSMark /* Assume that cells are originally numbered lexicographically */ 1110c569c6eSMark for (ck = 0; ck < pcells[2]; ++ck) { 1120c569c6eSMark for (cj = 0; cj < pcells[1]; ++cj) { 1130c569c6eSMark for (ci = 0; ci < pcells[0]; ++ci) { 1140c569c6eSMark 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; 1150c569c6eSMark 1160c569c6eSMark cellproc[offsets[pid]++] = cid; 1170c569c6eSMark } 1180c569c6eSMark } 1190c569c6eSMark } 1200c569c6eSMark } 1210c569c6eSMark } 1220c569c6eSMark } 1230c569c6eSMark } 1240c569c6eSMark } 1250c569c6eSMark } 1262c71b3e2SJacob Faibussowitsch for (np = 1; np < nparts; ++np) PetscCheckFalse(offsets[np] - offsets[np-1] != numVertices/nparts,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Offset %D != %D partition size", offsets[np], numVertices/nparts); 1279566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets)); 1289566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numVertices, cellproc, PETSC_OWN_POINTER, partition)); 1290c569c6eSMark PetscFunctionReturn(0); 1300c569c6eSMark } 1310c569c6eSMark 132abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 133abe9303eSLisandro Dalcin { 1340c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 135abe9303eSLisandro Dalcin MPI_Comm comm; 136abe9303eSLisandro Dalcin PetscInt np, *tpwgts = NULL, sumw = 0, numVerticesGlobal = 0; 137abe9303eSLisandro Dalcin PetscMPIInt size; 138abe9303eSLisandro Dalcin 139abe9303eSLisandro Dalcin PetscFunctionBegin; 1400c569c6eSMark if (p->useGrid) { 1419566063dSJacob Faibussowitsch PetscCall(PetscPartitionerPartition_Simple_Grid(part, nparts, numVertices, start, adjacency, vertSection, targetSection, partSection, partition)); 1420c569c6eSMark PetscFunctionReturn(0); 1430c569c6eSMark } 1449566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part,"PETSCPARTITIONERSIMPLE ignores vertex weights\n")); 1459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) part, &comm)); 1469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 147abe9303eSLisandro Dalcin if (targetSection) { 148*1c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm)); 1499566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts,&tpwgts)); 150abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection,np,&tpwgts[np])); 152abe9303eSLisandro Dalcin sumw += tpwgts[np]; 153abe9303eSLisandro Dalcin } 154d7cc930eSLisandro Dalcin if (sumw) { 155abe9303eSLisandro Dalcin PetscInt m,mp; 156abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np]*numVerticesGlobal)/sumw; 157abe9303eSLisandro Dalcin for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) { 158abe9303eSLisandro Dalcin if (m < tpwgts[np]) { m = tpwgts[np]; mp = np; } 159abe9303eSLisandro Dalcin sumw += tpwgts[np]; 160abe9303eSLisandro Dalcin } 161abe9303eSLisandro Dalcin if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw; 162abe9303eSLisandro Dalcin } 1639566063dSJacob Faibussowitsch if (!sumw) PetscCall(PetscFree(tpwgts)); 164abe9303eSLisandro Dalcin } 165abe9303eSLisandro Dalcin 1669566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition)); 167abe9303eSLisandro Dalcin if (size == 1) { 168abe9303eSLisandro Dalcin if (tpwgts) { 169abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1709566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, np, tpwgts[np])); 171abe9303eSLisandro Dalcin } 172abe9303eSLisandro Dalcin } else { 173abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1749566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np))); 175abe9303eSLisandro Dalcin } 176abe9303eSLisandro Dalcin } 177abe9303eSLisandro Dalcin } else { 178abe9303eSLisandro Dalcin if (tpwgts) { 179abe9303eSLisandro Dalcin Vec v; 180abe9303eSLisandro Dalcin PetscScalar *array; 181abe9303eSLisandro Dalcin PetscInt st,j; 182abe9303eSLisandro Dalcin PetscMPIInt rank; 183abe9303eSLisandro Dalcin 1849566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&v)); 1859566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v,numVertices,numVerticesGlobal)); 1869566063dSJacob Faibussowitsch PetscCall(VecSetType(v,VECSTANDARD)); 1879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 188abe9303eSLisandro Dalcin for (np = 0,st = 0; np < nparts; ++np) { 189abe9303eSLisandro Dalcin if (rank == np || (rank == size-1 && size < nparts && np >= size)) { 190abe9303eSLisandro Dalcin for (j = 0; j < tpwgts[np]; j++) { 1919566063dSJacob Faibussowitsch PetscCall(VecSetValue(v,st+j,np,INSERT_VALUES)); 192abe9303eSLisandro Dalcin } 193abe9303eSLisandro Dalcin } 194abe9303eSLisandro Dalcin st += tpwgts[np]; 195abe9303eSLisandro Dalcin } 1969566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(v)); 1979566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(v)); 1989566063dSJacob Faibussowitsch PetscCall(VecGetArray(v,&array)); 199abe9303eSLisandro Dalcin for (j = 0; j < numVertices; ++j) { 2009566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(partSection,PetscRealPart(array[j]),1)); 201abe9303eSLisandro Dalcin } 2029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&array)); 2039566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 204abe9303eSLisandro Dalcin } else { 205abe9303eSLisandro Dalcin PetscMPIInt rank; 206abe9303eSLisandro Dalcin PetscInt nvGlobal, *offsets, myFirst, myLast; 207abe9303eSLisandro Dalcin 2089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size+1,&offsets)); 209abe9303eSLisandro Dalcin offsets[0] = 0; 2109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm)); 211abe9303eSLisandro Dalcin for (np = 2; np <= size; np++) { 212abe9303eSLisandro Dalcin offsets[np] += offsets[np-1]; 213abe9303eSLisandro Dalcin } 214abe9303eSLisandro Dalcin nvGlobal = offsets[size]; 2159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 216abe9303eSLisandro Dalcin myFirst = offsets[rank]; 217abe9303eSLisandro Dalcin myLast = offsets[rank + 1] - 1; 2189566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets)); 219abe9303eSLisandro Dalcin if (numVertices) { 220abe9303eSLisandro Dalcin PetscInt firstPart = 0, firstLargePart = 0; 221abe9303eSLisandro Dalcin PetscInt lastPart = 0, lastLargePart = 0; 222abe9303eSLisandro Dalcin PetscInt rem = nvGlobal % nparts; 223abe9303eSLisandro Dalcin PetscInt pSmall = nvGlobal/nparts; 224abe9303eSLisandro Dalcin PetscInt pBig = nvGlobal/nparts + 1; 225abe9303eSLisandro Dalcin 226abe9303eSLisandro Dalcin if (rem) { 227abe9303eSLisandro Dalcin firstLargePart = myFirst / pBig; 228abe9303eSLisandro Dalcin lastLargePart = myLast / pBig; 229abe9303eSLisandro Dalcin 230abe9303eSLisandro Dalcin if (firstLargePart < rem) { 231abe9303eSLisandro Dalcin firstPart = firstLargePart; 232abe9303eSLisandro Dalcin } else { 233abe9303eSLisandro Dalcin firstPart = rem + (myFirst - (rem * pBig)) / pSmall; 234abe9303eSLisandro Dalcin } 235abe9303eSLisandro Dalcin if (lastLargePart < rem) { 236abe9303eSLisandro Dalcin lastPart = lastLargePart; 237abe9303eSLisandro Dalcin } else { 238abe9303eSLisandro Dalcin lastPart = rem + (myLast - (rem * pBig)) / pSmall; 239abe9303eSLisandro Dalcin } 240abe9303eSLisandro Dalcin } else { 241abe9303eSLisandro Dalcin firstPart = myFirst / (nvGlobal/nparts); 242abe9303eSLisandro Dalcin lastPart = myLast / (nvGlobal/nparts); 243abe9303eSLisandro Dalcin } 244abe9303eSLisandro Dalcin 245abe9303eSLisandro Dalcin for (np = firstPart; np <= lastPart; np++) { 246abe9303eSLisandro Dalcin PetscInt PartStart = np * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np); 247abe9303eSLisandro Dalcin PetscInt PartEnd = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1); 248abe9303eSLisandro Dalcin 249abe9303eSLisandro Dalcin PartStart = PetscMax(PartStart,myFirst); 250abe9303eSLisandro Dalcin PartEnd = PetscMin(PartEnd,myLast+1); 2519566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection,np,PartEnd-PartStart)); 252abe9303eSLisandro Dalcin } 253abe9303eSLisandro Dalcin } 254abe9303eSLisandro Dalcin } 255abe9303eSLisandro Dalcin } 2569566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts)); 257abe9303eSLisandro Dalcin PetscFunctionReturn(0); 258abe9303eSLisandro Dalcin } 259abe9303eSLisandro Dalcin 260abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part) 261abe9303eSLisandro Dalcin { 262abe9303eSLisandro Dalcin PetscFunctionBegin; 263abe9303eSLisandro Dalcin part->noGraph = PETSC_TRUE; 264abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_Simple; 2650c569c6eSMark part->ops->setfromoptions = PetscPartitionerSetFromOptions_Simple; 266abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_Simple; 267abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_Simple; 268abe9303eSLisandro Dalcin PetscFunctionReturn(0); 269abe9303eSLisandro Dalcin } 270abe9303eSLisandro Dalcin 271abe9303eSLisandro Dalcin /*MC 272abe9303eSLisandro Dalcin PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object 273abe9303eSLisandro Dalcin 274abe9303eSLisandro Dalcin Level: intermediate 275abe9303eSLisandro Dalcin 276abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 277abe9303eSLisandro Dalcin M*/ 278abe9303eSLisandro Dalcin 279abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part) 280abe9303eSLisandro Dalcin { 281abe9303eSLisandro Dalcin PetscPartitioner_Simple *p; 282abe9303eSLisandro Dalcin 283abe9303eSLisandro Dalcin PetscFunctionBegin; 284abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 2859566063dSJacob Faibussowitsch PetscCall(PetscNewLog(part, &p)); 2860c569c6eSMark p->gridDim = -1; 287abe9303eSLisandro Dalcin part->data = p; 288abe9303eSLisandro Dalcin 2899566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_Simple(part)); 290abe9303eSLisandro Dalcin PetscFunctionReturn(0); 291abe9303eSLisandro Dalcin } 292