xref: /petsc/src/dm/partitioner/impls/shell/partshell.c (revision abe9303e6325b68e0d8957680d58b261e7a295d5)
1*abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h>        /*I "petscpartitioner.h" I*/
2*abe9303eSLisandro Dalcin 
3*abe9303eSLisandro Dalcin typedef struct {
4*abe9303eSLisandro Dalcin   PetscSection section;   /* Sizes for each partition */
5*abe9303eSLisandro Dalcin   IS           partition; /* Points in each partition */
6*abe9303eSLisandro Dalcin   PetscBool    random;    /* Flag for a random partition */
7*abe9303eSLisandro Dalcin } PetscPartitioner_Shell;
8*abe9303eSLisandro Dalcin 
9*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerReset_Shell(PetscPartitioner part)
10*abe9303eSLisandro Dalcin {
11*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
12*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
13*abe9303eSLisandro Dalcin 
14*abe9303eSLisandro Dalcin   PetscFunctionBegin;
15*abe9303eSLisandro Dalcin   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
16*abe9303eSLisandro Dalcin   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
17*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
18*abe9303eSLisandro Dalcin }
19*abe9303eSLisandro Dalcin 
20*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
21*abe9303eSLisandro Dalcin {
22*abe9303eSLisandro Dalcin   PetscErrorCode ierr;
23*abe9303eSLisandro Dalcin 
24*abe9303eSLisandro Dalcin   PetscFunctionBegin;
25*abe9303eSLisandro Dalcin   ierr = PetscPartitionerReset_Shell(part);CHKERRQ(ierr);
26*abe9303eSLisandro Dalcin   ierr = PetscFree(part->data);CHKERRQ(ierr);
27*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
28*abe9303eSLisandro Dalcin }
29*abe9303eSLisandro Dalcin 
30*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Shell_ASCII(PetscPartitioner part, PetscViewer viewer)
31*abe9303eSLisandro Dalcin {
32*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
33*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
34*abe9303eSLisandro Dalcin 
35*abe9303eSLisandro Dalcin   PetscFunctionBegin;
36*abe9303eSLisandro Dalcin   if (p->random) {
37*abe9303eSLisandro Dalcin     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
38*abe9303eSLisandro Dalcin     ierr = PetscViewerASCIIPrintf(viewer, "using random partition\n");CHKERRQ(ierr);
39*abe9303eSLisandro Dalcin     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
40*abe9303eSLisandro Dalcin   }
41*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
42*abe9303eSLisandro Dalcin }
43*abe9303eSLisandro Dalcin 
44*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
45*abe9303eSLisandro Dalcin {
46*abe9303eSLisandro Dalcin   PetscBool      iascii;
47*abe9303eSLisandro Dalcin   PetscErrorCode ierr;
48*abe9303eSLisandro Dalcin 
49*abe9303eSLisandro Dalcin   PetscFunctionBegin;
50*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
51*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
52*abe9303eSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
53*abe9303eSLisandro Dalcin   if (iascii) {ierr = PetscPartitionerView_Shell_ASCII(part, viewer);CHKERRQ(ierr);}
54*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
55*abe9303eSLisandro Dalcin }
56*abe9303eSLisandro Dalcin 
57*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
58*abe9303eSLisandro Dalcin {
59*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
60*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
61*abe9303eSLisandro Dalcin 
62*abe9303eSLisandro Dalcin   PetscFunctionBegin;
63*abe9303eSLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitioner Shell Options");CHKERRQ(ierr);
64*abe9303eSLisandro Dalcin   ierr = PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &p->random, NULL);CHKERRQ(ierr);
65*abe9303eSLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
66*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
67*abe9303eSLisandro Dalcin }
68*abe9303eSLisandro Dalcin 
69*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
70*abe9303eSLisandro Dalcin {
71*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
72*abe9303eSLisandro Dalcin   PetscInt                np;
73*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
74*abe9303eSLisandro Dalcin 
75*abe9303eSLisandro Dalcin   PetscFunctionBegin;
76*abe9303eSLisandro Dalcin   if (p->random) {
77*abe9303eSLisandro Dalcin     PetscRandom r;
78*abe9303eSLisandro Dalcin     PetscInt   *sizes, *points, v, p;
79*abe9303eSLisandro Dalcin     PetscMPIInt rank;
80*abe9303eSLisandro Dalcin 
81*abe9303eSLisandro Dalcin     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) part), &rank);CHKERRQ(ierr);
82*abe9303eSLisandro Dalcin     ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
83*abe9303eSLisandro Dalcin     ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) nparts);CHKERRQ(ierr);
84*abe9303eSLisandro Dalcin     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
85*abe9303eSLisandro Dalcin     ierr = PetscCalloc2(nparts, &sizes, numVertices, &points);CHKERRQ(ierr);
86*abe9303eSLisandro Dalcin     for (v = 0; v < numVertices; ++v) {points[v] = v;}
87*abe9303eSLisandro Dalcin     for (p = 0; p < nparts; ++p) {sizes[p] = numVertices/nparts + (PetscInt) (p < numVertices % nparts);}
88*abe9303eSLisandro Dalcin     for (v = numVertices-1; v > 0; --v) {
89*abe9303eSLisandro Dalcin       PetscReal val;
90*abe9303eSLisandro Dalcin       PetscInt  w, tmp;
91*abe9303eSLisandro Dalcin 
92*abe9303eSLisandro Dalcin       ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) (v+1));CHKERRQ(ierr);
93*abe9303eSLisandro Dalcin       ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr);
94*abe9303eSLisandro Dalcin       w    = PetscFloorReal(val);
95*abe9303eSLisandro Dalcin       tmp       = points[v];
96*abe9303eSLisandro Dalcin       points[v] = points[w];
97*abe9303eSLisandro Dalcin       points[w] = tmp;
98*abe9303eSLisandro Dalcin     }
99*abe9303eSLisandro Dalcin     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
100*abe9303eSLisandro Dalcin     ierr = PetscPartitionerShellSetPartition(part, nparts, sizes, points);CHKERRQ(ierr);
101*abe9303eSLisandro Dalcin     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
102*abe9303eSLisandro Dalcin   }
103*abe9303eSLisandro Dalcin   if (!p->section) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()");
104*abe9303eSLisandro Dalcin   ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr);
105*abe9303eSLisandro Dalcin   if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
106*abe9303eSLisandro Dalcin   ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr);
107*abe9303eSLisandro Dalcin   if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
108*abe9303eSLisandro Dalcin   ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr);
109*abe9303eSLisandro Dalcin   *partition = p->partition;
110*abe9303eSLisandro Dalcin   ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr);
111*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
112*abe9303eSLisandro Dalcin }
113*abe9303eSLisandro Dalcin 
114*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
115*abe9303eSLisandro Dalcin {
116*abe9303eSLisandro Dalcin   PetscFunctionBegin;
117*abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE; /* PetscPartitionerShell cannot overload the partition call, so it is safe for now */
118*abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Shell;
119*abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
120*abe9303eSLisandro Dalcin   part->ops->reset          = PetscPartitionerReset_Shell;
121*abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Shell;
122*abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Shell;
123*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
124*abe9303eSLisandro Dalcin }
125*abe9303eSLisandro Dalcin 
126*abe9303eSLisandro Dalcin /*MC
127*abe9303eSLisandro Dalcin   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
128*abe9303eSLisandro Dalcin 
129*abe9303eSLisandro Dalcin   Level: intermediate
130*abe9303eSLisandro Dalcin 
131*abe9303eSLisandro Dalcin   Options Database Keys:
132*abe9303eSLisandro Dalcin .  -petscpartitioner_shell_random - Use a random partition
133*abe9303eSLisandro Dalcin 
134*abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
135*abe9303eSLisandro Dalcin M*/
136*abe9303eSLisandro Dalcin 
137*abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
138*abe9303eSLisandro Dalcin {
139*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p;
140*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
141*abe9303eSLisandro Dalcin 
142*abe9303eSLisandro Dalcin   PetscFunctionBegin;
143*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
144*abe9303eSLisandro Dalcin   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
145*abe9303eSLisandro Dalcin   part->data = p;
146*abe9303eSLisandro Dalcin 
147*abe9303eSLisandro Dalcin   ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr);
148*abe9303eSLisandro Dalcin   p->random = PETSC_FALSE;
149*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
150*abe9303eSLisandro Dalcin }
151*abe9303eSLisandro Dalcin 
152*abe9303eSLisandro Dalcin /*@C
153*abe9303eSLisandro Dalcin   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
154*abe9303eSLisandro Dalcin 
155*abe9303eSLisandro Dalcin   Collective on PetscPartitioner
156*abe9303eSLisandro Dalcin 
157*abe9303eSLisandro Dalcin   Input Parameters:
158*abe9303eSLisandro Dalcin + part   - The PetscPartitioner
159*abe9303eSLisandro Dalcin . size   - The number of partitions
160*abe9303eSLisandro Dalcin . sizes  - array of length size (or NULL) providing the number of points in each partition
161*abe9303eSLisandro Dalcin - points - array of length sum(sizes) (may be NULL iff sizes is NULL), a permutation of the points that groups those assigned to each partition in order (i.e., partition 0 first, partition 1 next, etc.)
162*abe9303eSLisandro Dalcin 
163*abe9303eSLisandro Dalcin   Level: developer
164*abe9303eSLisandro Dalcin 
165*abe9303eSLisandro Dalcin   Notes:
166*abe9303eSLisandro Dalcin     It is safe to free the sizes and points arrays after use in this routine.
167*abe9303eSLisandro Dalcin 
168*abe9303eSLisandro Dalcin .seealso DMPlexDistribute(), PetscPartitionerCreate()
169*abe9303eSLisandro Dalcin @*/
170*abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[])
171*abe9303eSLisandro Dalcin {
172*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
173*abe9303eSLisandro Dalcin   PetscInt                proc, numPoints;
174*abe9303eSLisandro Dalcin   PetscErrorCode          ierr;
175*abe9303eSLisandro Dalcin 
176*abe9303eSLisandro Dalcin   PetscFunctionBegin;
177*abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
178*abe9303eSLisandro Dalcin   if (sizes)  {PetscValidPointer(sizes, 3);}
179*abe9303eSLisandro Dalcin   if (points) {PetscValidPointer(points, 4);}
180*abe9303eSLisandro Dalcin   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
181*abe9303eSLisandro Dalcin   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
182*abe9303eSLisandro Dalcin   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr);
183*abe9303eSLisandro Dalcin   ierr = PetscSectionSetChart(p->section, 0, size);CHKERRQ(ierr);
184*abe9303eSLisandro Dalcin   if (sizes) {
185*abe9303eSLisandro Dalcin     for (proc = 0; proc < size; ++proc) {
186*abe9303eSLisandro Dalcin       ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr);
187*abe9303eSLisandro Dalcin     }
188*abe9303eSLisandro Dalcin   }
189*abe9303eSLisandro Dalcin   ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr);
190*abe9303eSLisandro Dalcin   ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr);
191*abe9303eSLisandro Dalcin   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr);
192*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
193*abe9303eSLisandro Dalcin }
194*abe9303eSLisandro Dalcin 
195*abe9303eSLisandro Dalcin /*@
196*abe9303eSLisandro Dalcin   PetscPartitionerShellSetRandom - Set the flag to use a random partition
197*abe9303eSLisandro Dalcin 
198*abe9303eSLisandro Dalcin   Collective on PetscPartitioner
199*abe9303eSLisandro Dalcin 
200*abe9303eSLisandro Dalcin   Input Parameters:
201*abe9303eSLisandro Dalcin + part   - The PetscPartitioner
202*abe9303eSLisandro Dalcin - random - The flag to use a random partition
203*abe9303eSLisandro Dalcin 
204*abe9303eSLisandro Dalcin   Level: intermediate
205*abe9303eSLisandro Dalcin 
206*abe9303eSLisandro Dalcin .seealso PetscPartitionerShellGetRandom(), PetscPartitionerCreate()
207*abe9303eSLisandro Dalcin @*/
208*abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random)
209*abe9303eSLisandro Dalcin {
210*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
211*abe9303eSLisandro Dalcin 
212*abe9303eSLisandro Dalcin   PetscFunctionBegin;
213*abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
214*abe9303eSLisandro Dalcin   p->random = random;
215*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
216*abe9303eSLisandro Dalcin }
217*abe9303eSLisandro Dalcin 
218*abe9303eSLisandro Dalcin /*@
219*abe9303eSLisandro Dalcin   PetscPartitionerShellGetRandom - get the flag to use a random partition
220*abe9303eSLisandro Dalcin 
221*abe9303eSLisandro Dalcin   Collective on PetscPartitioner
222*abe9303eSLisandro Dalcin 
223*abe9303eSLisandro Dalcin   Input Parameter:
224*abe9303eSLisandro Dalcin . part   - The PetscPartitioner
225*abe9303eSLisandro Dalcin 
226*abe9303eSLisandro Dalcin   Output Parameter:
227*abe9303eSLisandro Dalcin . random - The flag to use a random partition
228*abe9303eSLisandro Dalcin 
229*abe9303eSLisandro Dalcin   Level: intermediate
230*abe9303eSLisandro Dalcin 
231*abe9303eSLisandro Dalcin .seealso PetscPartitionerShellSetRandom(), PetscPartitionerCreate()
232*abe9303eSLisandro Dalcin @*/
233*abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random)
234*abe9303eSLisandro Dalcin {
235*abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
236*abe9303eSLisandro Dalcin 
237*abe9303eSLisandro Dalcin   PetscFunctionBegin;
238*abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
239*abe9303eSLisandro Dalcin   PetscValidPointer(random, 2);
240*abe9303eSLisandro Dalcin   *random = p->random;
241*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
242*abe9303eSLisandro Dalcin }
243