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