1abe9303eSLisandro Dalcin #include <petscvec.h> 2abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/ 3abe9303eSLisandro Dalcin 4abe9303eSLisandro Dalcin typedef struct { 5abe9303eSLisandro Dalcin PetscInt dummy; 6abe9303eSLisandro Dalcin } PetscPartitioner_Simple; 7abe9303eSLisandro Dalcin 8abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part) 9abe9303eSLisandro Dalcin { 10abe9303eSLisandro Dalcin PetscErrorCode ierr; 11abe9303eSLisandro Dalcin 12abe9303eSLisandro Dalcin PetscFunctionBegin; 13abe9303eSLisandro Dalcin ierr = PetscFree(part->data);CHKERRQ(ierr); 14abe9303eSLisandro Dalcin PetscFunctionReturn(0); 15abe9303eSLisandro Dalcin } 16abe9303eSLisandro Dalcin 17abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple_ASCII(PetscPartitioner part, PetscViewer viewer) 18abe9303eSLisandro Dalcin { 19abe9303eSLisandro Dalcin PetscFunctionBegin; 20abe9303eSLisandro Dalcin PetscFunctionReturn(0); 21abe9303eSLisandro Dalcin } 22abe9303eSLisandro Dalcin 23abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer) 24abe9303eSLisandro Dalcin { 25abe9303eSLisandro Dalcin PetscBool iascii; 26abe9303eSLisandro Dalcin PetscErrorCode ierr; 27abe9303eSLisandro Dalcin 28abe9303eSLisandro Dalcin PetscFunctionBegin; 29abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 30abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 31abe9303eSLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 32abe9303eSLisandro Dalcin if (iascii) {ierr = PetscPartitionerView_Simple_ASCII(part, viewer);CHKERRQ(ierr);} 33abe9303eSLisandro Dalcin PetscFunctionReturn(0); 34abe9303eSLisandro Dalcin } 35abe9303eSLisandro Dalcin 36abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 37abe9303eSLisandro Dalcin { 38abe9303eSLisandro Dalcin MPI_Comm comm; 39abe9303eSLisandro Dalcin PetscInt np, *tpwgts = NULL, sumw = 0, numVerticesGlobal = 0; 40abe9303eSLisandro Dalcin PetscMPIInt size; 41abe9303eSLisandro Dalcin PetscErrorCode ierr; 42abe9303eSLisandro Dalcin 43abe9303eSLisandro Dalcin PetscFunctionBegin; 44abe9303eSLisandro Dalcin if (vertSection) { ierr = PetscInfo(part,"PETSCPARTITIONERSIMPLE ignores vertex weights\n");CHKERRQ(ierr); } 45abe9303eSLisandro Dalcin comm = PetscObjectComm((PetscObject)part); 46abe9303eSLisandro Dalcin ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 47abe9303eSLisandro Dalcin if (targetSection) { 48abe9303eSLisandro Dalcin ierr = MPIU_Allreduce(&numVertices, &numVerticesGlobal, 1, MPIU_INT, MPI_SUM, comm);CHKERRQ(ierr); 49abe9303eSLisandro Dalcin ierr = PetscCalloc1(nparts,&tpwgts);CHKERRQ(ierr); 50abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 51abe9303eSLisandro Dalcin ierr = PetscSectionGetDof(targetSection,np,&tpwgts[np]);CHKERRQ(ierr); 52abe9303eSLisandro Dalcin sumw += tpwgts[np]; 53abe9303eSLisandro Dalcin } 54*d7cc930eSLisandro Dalcin if (sumw) { 55abe9303eSLisandro Dalcin PetscInt m,mp; 56abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) tpwgts[np] = (tpwgts[np]*numVerticesGlobal)/sumw; 57abe9303eSLisandro Dalcin for (np = 0, m = -1, mp = 0, sumw = 0; np < nparts; ++np) { 58abe9303eSLisandro Dalcin if (m < tpwgts[np]) { m = tpwgts[np]; mp = np; } 59abe9303eSLisandro Dalcin sumw += tpwgts[np]; 60abe9303eSLisandro Dalcin } 61abe9303eSLisandro Dalcin if (sumw != numVerticesGlobal) tpwgts[mp] += numVerticesGlobal - sumw; 62abe9303eSLisandro Dalcin } 63*d7cc930eSLisandro Dalcin if (!sumw) {ierr = PetscFree(tpwgts);CHKERRQ(ierr);} 64abe9303eSLisandro Dalcin } 65abe9303eSLisandro Dalcin 66abe9303eSLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr); 67abe9303eSLisandro Dalcin if (size == 1) { 68abe9303eSLisandro Dalcin if (tpwgts) { 69abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 70abe9303eSLisandro Dalcin ierr = PetscSectionSetDof(partSection, np, tpwgts[np]);CHKERRQ(ierr); 71abe9303eSLisandro Dalcin } 72abe9303eSLisandro Dalcin } else { 73abe9303eSLisandro Dalcin for (np = 0; np < nparts; ++np) { 74abe9303eSLisandro Dalcin ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr); 75abe9303eSLisandro Dalcin } 76abe9303eSLisandro Dalcin } 77abe9303eSLisandro Dalcin } else { 78abe9303eSLisandro Dalcin if (tpwgts) { 79abe9303eSLisandro Dalcin Vec v; 80abe9303eSLisandro Dalcin PetscScalar *array; 81abe9303eSLisandro Dalcin PetscInt st,j; 82abe9303eSLisandro Dalcin PetscMPIInt rank; 83abe9303eSLisandro Dalcin 84abe9303eSLisandro Dalcin ierr = VecCreate(comm,&v);CHKERRQ(ierr); 85abe9303eSLisandro Dalcin ierr = VecSetSizes(v,numVertices,numVerticesGlobal);CHKERRQ(ierr); 86abe9303eSLisandro Dalcin ierr = VecSetType(v,VECSTANDARD);CHKERRQ(ierr); 87abe9303eSLisandro Dalcin ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 88abe9303eSLisandro Dalcin for (np = 0,st = 0; np < nparts; ++np) { 89abe9303eSLisandro Dalcin if (rank == np || (rank == size-1 && size < nparts && np >= size)) { 90abe9303eSLisandro Dalcin for (j = 0; j < tpwgts[np]; j++) { 91abe9303eSLisandro Dalcin ierr = VecSetValue(v,st+j,np,INSERT_VALUES);CHKERRQ(ierr); 92abe9303eSLisandro Dalcin } 93abe9303eSLisandro Dalcin } 94abe9303eSLisandro Dalcin st += tpwgts[np]; 95abe9303eSLisandro Dalcin } 96abe9303eSLisandro Dalcin ierr = VecAssemblyBegin(v);CHKERRQ(ierr); 97abe9303eSLisandro Dalcin ierr = VecAssemblyEnd(v);CHKERRQ(ierr); 98abe9303eSLisandro Dalcin ierr = VecGetArray(v,&array);CHKERRQ(ierr); 99abe9303eSLisandro Dalcin for (j = 0; j < numVertices; ++j) { 100abe9303eSLisandro Dalcin ierr = PetscSectionAddDof(partSection,PetscRealPart(array[j]),1);CHKERRQ(ierr); 101abe9303eSLisandro Dalcin } 102abe9303eSLisandro Dalcin ierr = VecRestoreArray(v,&array);CHKERRQ(ierr); 103abe9303eSLisandro Dalcin ierr = VecDestroy(&v);CHKERRQ(ierr); 104abe9303eSLisandro Dalcin } else { 105abe9303eSLisandro Dalcin PetscMPIInt rank; 106abe9303eSLisandro Dalcin PetscInt nvGlobal, *offsets, myFirst, myLast; 107abe9303eSLisandro Dalcin 108abe9303eSLisandro Dalcin ierr = PetscMalloc1(size+1,&offsets);CHKERRQ(ierr); 109abe9303eSLisandro Dalcin offsets[0] = 0; 110abe9303eSLisandro Dalcin ierr = MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm);CHKERRQ(ierr); 111abe9303eSLisandro Dalcin for (np = 2; np <= size; np++) { 112abe9303eSLisandro Dalcin offsets[np] += offsets[np-1]; 113abe9303eSLisandro Dalcin } 114abe9303eSLisandro Dalcin nvGlobal = offsets[size]; 115abe9303eSLisandro Dalcin ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 116abe9303eSLisandro Dalcin myFirst = offsets[rank]; 117abe9303eSLisandro Dalcin myLast = offsets[rank + 1] - 1; 118abe9303eSLisandro Dalcin ierr = PetscFree(offsets);CHKERRQ(ierr); 119abe9303eSLisandro Dalcin if (numVertices) { 120abe9303eSLisandro Dalcin PetscInt firstPart = 0, firstLargePart = 0; 121abe9303eSLisandro Dalcin PetscInt lastPart = 0, lastLargePart = 0; 122abe9303eSLisandro Dalcin PetscInt rem = nvGlobal % nparts; 123abe9303eSLisandro Dalcin PetscInt pSmall = nvGlobal/nparts; 124abe9303eSLisandro Dalcin PetscInt pBig = nvGlobal/nparts + 1; 125abe9303eSLisandro Dalcin 126abe9303eSLisandro Dalcin if (rem) { 127abe9303eSLisandro Dalcin firstLargePart = myFirst / pBig; 128abe9303eSLisandro Dalcin lastLargePart = myLast / pBig; 129abe9303eSLisandro Dalcin 130abe9303eSLisandro Dalcin if (firstLargePart < rem) { 131abe9303eSLisandro Dalcin firstPart = firstLargePart; 132abe9303eSLisandro Dalcin } else { 133abe9303eSLisandro Dalcin firstPart = rem + (myFirst - (rem * pBig)) / pSmall; 134abe9303eSLisandro Dalcin } 135abe9303eSLisandro Dalcin if (lastLargePart < rem) { 136abe9303eSLisandro Dalcin lastPart = lastLargePart; 137abe9303eSLisandro Dalcin } else { 138abe9303eSLisandro Dalcin lastPart = rem + (myLast - (rem * pBig)) / pSmall; 139abe9303eSLisandro Dalcin } 140abe9303eSLisandro Dalcin } else { 141abe9303eSLisandro Dalcin firstPart = myFirst / (nvGlobal/nparts); 142abe9303eSLisandro Dalcin lastPart = myLast / (nvGlobal/nparts); 143abe9303eSLisandro Dalcin } 144abe9303eSLisandro Dalcin 145abe9303eSLisandro Dalcin for (np = firstPart; np <= lastPart; np++) { 146abe9303eSLisandro Dalcin PetscInt PartStart = np * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np); 147abe9303eSLisandro Dalcin PetscInt PartEnd = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1); 148abe9303eSLisandro Dalcin 149abe9303eSLisandro Dalcin PartStart = PetscMax(PartStart,myFirst); 150abe9303eSLisandro Dalcin PartEnd = PetscMin(PartEnd,myLast+1); 151abe9303eSLisandro Dalcin ierr = PetscSectionSetDof(partSection,np,PartEnd-PartStart);CHKERRQ(ierr); 152abe9303eSLisandro Dalcin } 153abe9303eSLisandro Dalcin } 154abe9303eSLisandro Dalcin } 155abe9303eSLisandro Dalcin } 156abe9303eSLisandro Dalcin ierr = PetscFree(tpwgts);CHKERRQ(ierr); 157abe9303eSLisandro Dalcin PetscFunctionReturn(0); 158abe9303eSLisandro Dalcin } 159abe9303eSLisandro Dalcin 160abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part) 161abe9303eSLisandro Dalcin { 162abe9303eSLisandro Dalcin PetscFunctionBegin; 163abe9303eSLisandro Dalcin part->noGraph = PETSC_TRUE; 164abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_Simple; 165abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_Simple; 166abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_Simple; 167abe9303eSLisandro Dalcin PetscFunctionReturn(0); 168abe9303eSLisandro Dalcin } 169abe9303eSLisandro Dalcin 170abe9303eSLisandro Dalcin /*MC 171abe9303eSLisandro Dalcin PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object 172abe9303eSLisandro Dalcin 173abe9303eSLisandro Dalcin Level: intermediate 174abe9303eSLisandro Dalcin 175abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType() 176abe9303eSLisandro Dalcin M*/ 177abe9303eSLisandro Dalcin 178abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part) 179abe9303eSLisandro Dalcin { 180abe9303eSLisandro Dalcin PetscPartitioner_Simple *p; 181abe9303eSLisandro Dalcin PetscErrorCode ierr; 182abe9303eSLisandro Dalcin 183abe9303eSLisandro Dalcin PetscFunctionBegin; 184abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 185abe9303eSLisandro Dalcin ierr = PetscNewLog(part, &p);CHKERRQ(ierr); 186abe9303eSLisandro Dalcin part->data = p; 187abe9303eSLisandro Dalcin 188abe9303eSLisandro Dalcin ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr); 189abe9303eSLisandro Dalcin PetscFunctionReturn(0); 190abe9303eSLisandro Dalcin } 191