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 12abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part) 13abe9303eSLisandro Dalcin { 14abe9303eSLisandro Dalcin PetscFunctionBegin; 159566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 16abe9303eSLisandro Dalcin PetscFunctionReturn(0); 17abe9303eSLisandro Dalcin } 18abe9303eSLisandro Dalcin 19abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple_ASCII(PetscPartitioner part, PetscViewer viewer) 20abe9303eSLisandro Dalcin { 21abe9303eSLisandro Dalcin PetscFunctionBegin; 22abe9303eSLisandro Dalcin PetscFunctionReturn(0); 23abe9303eSLisandro Dalcin } 24abe9303eSLisandro Dalcin 25abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer) 26abe9303eSLisandro Dalcin { 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)); 34abe9303eSLisandro Dalcin PetscFunctionReturn(0); 35abe9303eSLisandro Dalcin } 36abe9303eSLisandro Dalcin 370c569c6eSMark static PetscErrorCode PetscPartitionerSetFromOptions_Simple(PetscOptionItems *PetscOptionsObject, PetscPartitioner part) 380c569c6eSMark { 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)); 480c569c6eSMark if (flg) {p->useGrid = PETSC_TRUE; p->gridDim = num;} 490c569c6eSMark num = 3; 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-petscpartitioner_simple_process_grid", "Number of local processes in each dimension for a given node", "", p->processGrid, &num, &flg)); 510c569c6eSMark if (flg) { 520c569c6eSMark p->useGrid = PETSC_TRUE; 530c569c6eSMark if (p->gridDim < 0) p->gridDim = num; 5463a3b9bcSJacob 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); 550c569c6eSMark } 56d0609cedSBarry Smith PetscOptionsHeadEnd(); 570c569c6eSMark PetscFunctionReturn(0); 580c569c6eSMark } 590c569c6eSMark 600c569c6eSMark static PetscErrorCode PetscPartitionerPartition_Simple_Grid(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 610c569c6eSMark { 620c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 630c569c6eSMark const PetscInt *nodes = p->nodeGrid; 640c569c6eSMark const PetscInt *procs = p->processGrid; 650c569c6eSMark PetscInt *cellproc, *offsets, cells[3] = {1, 1, 1}, pcells[3] = {1, 1, 1}; 66bd026e97SJed Brown PetscInt Np = 1, Nr, np, nk, nj, ni, pk, pj, pi, ck, cj, ci, i; 670c569c6eSMark MPI_Comm comm; 680c569c6eSMark PetscMPIInt size; 690c569c6eSMark 700c569c6eSMark PetscFunctionBegin; 719566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores vertex weights when using grid partition\n")); 729566063dSJacob Faibussowitsch if (targetSection) PetscCall(PetscInfo(part, "PETSCPARTITIONERSIMPLE ignores partition weights when using grid partition\n")); 739566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) part, &comm)); 749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 750c569c6eSMark /* Check grid */ 760c569c6eSMark for (i = 0; i < 3; ++i) Np *= nodes[i]*procs[i]; 7763a3b9bcSJacob Faibussowitsch PetscCheck(nparts == Np,comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %" PetscInt_FMT " grid size", nparts, Np); 7863a3b9bcSJacob Faibussowitsch PetscCheck(nparts == size,comm, PETSC_ERR_ARG_INCOMP, "Number of partitions %" PetscInt_FMT " != %d processes", nparts, size); 7963a3b9bcSJacob 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); 800c569c6eSMark for (i = 0; i < p->gridDim; ++i) cells[i] = nodes[i]*procs[i]; 810c569c6eSMark Nr = numVertices / nparts; 820c569c6eSMark while (Nr > 1) { 830c569c6eSMark for (i = 0; i < p->gridDim; ++i) { 840c569c6eSMark cells[i] *= 2; 850c569c6eSMark Nr /= 2; 860c569c6eSMark } 870c569c6eSMark } 88*1dca8a05SBarry Smith PetscCheck(!numVertices || Nr == 1,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Odd number of cells %" PetscInt_FMT ". Must be nprocs*2^k", numVertices); 890c569c6eSMark for (i = 0; i < p->gridDim; ++i) { 9063a3b9bcSJacob 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]); 91bd026e97SJed Brown pcells[i] = cells[i] / (nodes[i]*procs[i]); 920c569c6eSMark } 930c569c6eSMark /* Compute sizes */ 949566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, numVertices/nparts)); 959566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(partSection)); 969566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts, &offsets)); 979566063dSJacob Faibussowitsch for (np = 0; np < nparts; ++np) PetscCall(PetscSectionGetOffset(partSection, np, &offsets[np])); 980c569c6eSMark if (!numVertices) pcells[0] = pcells[1] = pcells[2] = 0; 990c569c6eSMark /* Compute partition */ 1009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numVertices, &cellproc)); 1010c569c6eSMark for (nk = 0; nk < nodes[2]; ++nk) { 1020c569c6eSMark for (nj = 0; nj < nodes[1]; ++nj) { 1030c569c6eSMark for (ni = 0; ni < nodes[0]; ++ni) { 1040c569c6eSMark const PetscInt nid = (nk*nodes[1] + nj)*nodes[0] + ni; 1050c569c6eSMark 1060c569c6eSMark for (pk = 0; pk < procs[2]; ++pk) { 1070c569c6eSMark for (pj = 0; pj < procs[1]; ++pj) { 1080c569c6eSMark for (pi = 0; pi < procs[0]; ++pi) { 1090c569c6eSMark const PetscInt pid = ((nid*procs[2] + pk)*procs[1] + pj)*procs[0] + pi; 1100c569c6eSMark 1110c569c6eSMark /* Assume that cells are originally numbered lexicographically */ 1120c569c6eSMark for (ck = 0; ck < pcells[2]; ++ck) { 1130c569c6eSMark for (cj = 0; cj < pcells[1]; ++cj) { 1140c569c6eSMark for (ci = 0; ci < pcells[0]; ++ci) { 1150c569c6eSMark 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; 1160c569c6eSMark 1170c569c6eSMark cellproc[offsets[pid]++] = cid; 1180c569c6eSMark } 1190c569c6eSMark } 1200c569c6eSMark } 1210c569c6eSMark } 1220c569c6eSMark } 1230c569c6eSMark } 1240c569c6eSMark } 1250c569c6eSMark } 1260c569c6eSMark } 127*1dca8a05SBarry 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); 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets)); 1299566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numVertices, cellproc, PETSC_OWN_POINTER, partition)); 1300c569c6eSMark PetscFunctionReturn(0); 1310c569c6eSMark } 1320c569c6eSMark 133abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 134abe9303eSLisandro Dalcin { 1350c569c6eSMark PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data; 136abe9303eSLisandro Dalcin MPI_Comm comm; 137abe9303eSLisandro Dalcin PetscInt np, *tpwgts = NULL, sumw = 0, numVerticesGlobal = 0; 138abe9303eSLisandro Dalcin PetscMPIInt size; 139abe9303eSLisandro Dalcin 140abe9303eSLisandro Dalcin PetscFunctionBegin; 1410c569c6eSMark if (p->useGrid) { 1429566063dSJacob Faibussowitsch PetscCall(PetscPartitionerPartition_Simple_Grid(part, nparts, numVertices, start, adjacency, vertSection, targetSection, partSection, partition)); 1430c569c6eSMark PetscFunctionReturn(0); 1440c569c6eSMark } 1459566063dSJacob Faibussowitsch if (vertSection) PetscCall(PetscInfo(part,"PETSCPARTITIONERSIMPLE ignores vertex weights\n")); 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) part, &comm)); 1479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 148abe9303eSLisandro Dalcin if (targetSection) { 1491c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm)); 1509566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(nparts,&tpwgts)); 151abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(targetSection,np,&tpwgts[np])); 153abe9303eSLisandro Dalcin sumw += tpwgts[np]; 154abe9303eSLisandro Dalcin } 155d7cc930eSLisandro Dalcin if (sumw) { 156abe9303eSLisandro Dalcin PetscInt m,mp; 157abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np]*numVerticesGlobal)/sumw; 158abe9303eSLisandro Dalcin for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) { 159abe9303eSLisandro Dalcin if (m < tpwgts[np]) { m = tpwgts[np]; mp = np; } 160abe9303eSLisandro Dalcin sumw += tpwgts[np]; 161abe9303eSLisandro Dalcin } 162abe9303eSLisandro Dalcin if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw; 163abe9303eSLisandro Dalcin } 1649566063dSJacob Faibussowitsch if (!sumw) PetscCall(PetscFree(tpwgts)); 165abe9303eSLisandro Dalcin } 166abe9303eSLisandro Dalcin 1679566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition)); 168abe9303eSLisandro Dalcin if (size == 1) { 169abe9303eSLisandro Dalcin if (tpwgts) { 170abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1719566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, np, tpwgts[np])); 172abe9303eSLisandro Dalcin } 173abe9303eSLisandro Dalcin } else { 174abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 1759566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np))); 176abe9303eSLisandro Dalcin } 177abe9303eSLisandro Dalcin } 178abe9303eSLisandro Dalcin } else { 179abe9303eSLisandro Dalcin if (tpwgts) { 180abe9303eSLisandro Dalcin Vec v; 181abe9303eSLisandro Dalcin PetscScalar *array; 182abe9303eSLisandro Dalcin PetscInt st,j; 183abe9303eSLisandro Dalcin PetscMPIInt rank; 184abe9303eSLisandro Dalcin 1859566063dSJacob Faibussowitsch PetscCall(VecCreate(comm,&v)); 1869566063dSJacob Faibussowitsch PetscCall(VecSetSizes(v,numVertices,numVerticesGlobal)); 1879566063dSJacob Faibussowitsch PetscCall(VecSetType(v,VECSTANDARD)); 1889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 189abe9303eSLisandro Dalcin for (np = 0,st = 0; np < nparts; ++np) { 190abe9303eSLisandro Dalcin if (rank == np || (rank == size-1 && size < nparts && np >= size)) { 191abe9303eSLisandro Dalcin for (j = 0; j < tpwgts[np]; j++) { 1929566063dSJacob Faibussowitsch PetscCall(VecSetValue(v,st+j,np,INSERT_VALUES)); 193abe9303eSLisandro Dalcin } 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)); 200abe9303eSLisandro Dalcin for (j = 0; j < numVertices; ++j) { 2019566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(partSection,PetscRealPart(array[j]),1)); 202abe9303eSLisandro Dalcin } 2039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v,&array)); 2049566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 205abe9303eSLisandro Dalcin } else { 206abe9303eSLisandro Dalcin PetscMPIInt rank; 207abe9303eSLisandro Dalcin PetscInt nvGlobal, *offsets, myFirst, myLast; 208abe9303eSLisandro Dalcin 2099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size+1,&offsets)); 210abe9303eSLisandro Dalcin offsets[0] = 0; 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm)); 212abe9303eSLisandro Dalcin for (np = 2; np <= size; np++) { 213abe9303eSLisandro Dalcin offsets[np] += offsets[np-1]; 214abe9303eSLisandro Dalcin } 215abe9303eSLisandro Dalcin nvGlobal = offsets[size]; 2169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 217abe9303eSLisandro Dalcin myFirst = offsets[rank]; 218abe9303eSLisandro Dalcin myLast = offsets[rank + 1] - 1; 2199566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets)); 220abe9303eSLisandro Dalcin if (numVertices) { 221abe9303eSLisandro Dalcin PetscInt firstPart = 0, firstLargePart = 0; 222abe9303eSLisandro Dalcin PetscInt lastPart = 0, lastLargePart = 0; 223abe9303eSLisandro Dalcin PetscInt rem = nvGlobal % nparts; 224abe9303eSLisandro Dalcin PetscInt pSmall = nvGlobal/nparts; 225abe9303eSLisandro Dalcin PetscInt pBig = nvGlobal/nparts + 1; 226abe9303eSLisandro Dalcin 227abe9303eSLisandro Dalcin if (rem) { 228abe9303eSLisandro Dalcin firstLargePart = myFirst / pBig; 229abe9303eSLisandro Dalcin lastLargePart = myLast / pBig; 230abe9303eSLisandro Dalcin 231abe9303eSLisandro Dalcin if (firstLargePart < rem) { 232abe9303eSLisandro Dalcin firstPart = firstLargePart; 233abe9303eSLisandro Dalcin } else { 234abe9303eSLisandro Dalcin firstPart = rem + (myFirst - (rem * pBig)) / pSmall; 235abe9303eSLisandro Dalcin } 236abe9303eSLisandro Dalcin if (lastLargePart < rem) { 237abe9303eSLisandro Dalcin lastPart = lastLargePart; 238abe9303eSLisandro Dalcin } else { 239abe9303eSLisandro Dalcin lastPart = rem + (myLast - (rem * pBig)) / pSmall; 240abe9303eSLisandro Dalcin } 241abe9303eSLisandro Dalcin } else { 242abe9303eSLisandro Dalcin firstPart = myFirst / (nvGlobal/nparts); 243abe9303eSLisandro Dalcin lastPart = myLast / (nvGlobal/nparts); 244abe9303eSLisandro Dalcin } 245abe9303eSLisandro Dalcin 246abe9303eSLisandro Dalcin for (np = firstPart; np <= lastPart; np++) { 247abe9303eSLisandro Dalcin PetscInt PartStart = np * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np); 248abe9303eSLisandro Dalcin PetscInt PartEnd = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1); 249abe9303eSLisandro Dalcin 250abe9303eSLisandro Dalcin PartStart = PetscMax(PartStart,myFirst); 251abe9303eSLisandro Dalcin PartEnd = PetscMin(PartEnd,myLast+1); 2529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection,np,PartEnd-PartStart)); 253abe9303eSLisandro Dalcin } 254abe9303eSLisandro Dalcin } 255abe9303eSLisandro Dalcin } 256abe9303eSLisandro Dalcin } 2579566063dSJacob Faibussowitsch PetscCall(PetscFree(tpwgts)); 258abe9303eSLisandro Dalcin PetscFunctionReturn(0); 259abe9303eSLisandro Dalcin } 260abe9303eSLisandro Dalcin 261abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part) 262abe9303eSLisandro Dalcin { 263abe9303eSLisandro Dalcin PetscFunctionBegin; 264abe9303eSLisandro Dalcin part->noGraph = PETSC_TRUE; 265abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_Simple; 2660c569c6eSMark part->ops->setfromoptions = PetscPartitionerSetFromOptions_Simple; 267abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_Simple; 268abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_Simple; 269abe9303eSLisandro Dalcin PetscFunctionReturn(0); 270abe9303eSLisandro Dalcin } 271abe9303eSLisandro Dalcin 272abe9303eSLisandro Dalcin /*MC 273abe9303eSLisandro Dalcin PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object 274abe9303eSLisandro Dalcin 275abe9303eSLisandro Dalcin Level: intermediate 276abe9303eSLisandro Dalcin 277abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 278abe9303eSLisandro Dalcin M*/ 279abe9303eSLisandro Dalcin 280abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part) 281abe9303eSLisandro Dalcin { 282abe9303eSLisandro Dalcin PetscPartitioner_Simple *p; 283abe9303eSLisandro Dalcin 284abe9303eSLisandro Dalcin PetscFunctionBegin; 285abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 2869566063dSJacob Faibussowitsch PetscCall(PetscNewLog(part, &p)); 2870c569c6eSMark p->gridDim = -1; 288abe9303eSLisandro Dalcin part->data = p; 289abe9303eSLisandro Dalcin 2909566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_Simple(part)); 291abe9303eSLisandro Dalcin PetscFunctionReturn(0); 292abe9303eSLisandro Dalcin } 293