xref: /petsc/src/dm/partitioner/impls/shell/partshell.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h>        /*I "petscpartitioner.h" I*/
2abe9303eSLisandro Dalcin 
3abe9303eSLisandro Dalcin typedef struct {
4abe9303eSLisandro Dalcin   PetscSection section;   /* Sizes for each partition */
5abe9303eSLisandro Dalcin   IS           partition; /* Points in each partition */
6abe9303eSLisandro Dalcin   PetscBool    random;    /* Flag for a random partition */
7abe9303eSLisandro Dalcin } PetscPartitioner_Shell;
8abe9303eSLisandro Dalcin 
9abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerReset_Shell(PetscPartitioner part)
10abe9303eSLisandro Dalcin {
11abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
12abe9303eSLisandro Dalcin 
13abe9303eSLisandro Dalcin   PetscFunctionBegin;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&p->section));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&p->partition));
16abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
17abe9303eSLisandro Dalcin }
18abe9303eSLisandro Dalcin 
19abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
20abe9303eSLisandro Dalcin {
21abe9303eSLisandro Dalcin   PetscFunctionBegin;
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerReset_Shell(part));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(part->data));
24abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
25abe9303eSLisandro Dalcin }
26abe9303eSLisandro Dalcin 
27abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Shell_ASCII(PetscPartitioner part, PetscViewer viewer)
28abe9303eSLisandro Dalcin {
29abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
30abe9303eSLisandro Dalcin 
31abe9303eSLisandro Dalcin   PetscFunctionBegin;
32abe9303eSLisandro Dalcin   if (p->random) {
335f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPushTab(viewer));
345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPrintf(viewer, "using random partition\n"));
355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerASCIIPopTab(viewer));
36abe9303eSLisandro Dalcin   }
37abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
38abe9303eSLisandro Dalcin }
39abe9303eSLisandro Dalcin 
40abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
41abe9303eSLisandro Dalcin {
42abe9303eSLisandro Dalcin   PetscBool      iascii;
43abe9303eSLisandro Dalcin 
44abe9303eSLisandro Dalcin   PetscFunctionBegin;
45abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
46abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii));
485f80ce2aSJacob Faibussowitsch   if (iascii) CHKERRQ(PetscPartitionerView_Shell_ASCII(part, viewer));
49abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
50abe9303eSLisandro Dalcin }
51abe9303eSLisandro Dalcin 
52abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
53abe9303eSLisandro Dalcin {
54d7cc930eSLisandro Dalcin   PetscBool      random = PETSC_FALSE, set;
55abe9303eSLisandro Dalcin 
56abe9303eSLisandro Dalcin   PetscFunctionBegin;
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHead(PetscOptionsObject, "PetscPartitioner Shell Options"));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerShellGetRandom(part, &random));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &random, &set));
605f80ce2aSJacob Faibussowitsch   if (set) CHKERRQ(PetscPartitionerShellSetRandom(part, random));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsTail());
62abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
63abe9303eSLisandro Dalcin }
64abe9303eSLisandro Dalcin 
65abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
66abe9303eSLisandro Dalcin {
67abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
68abe9303eSLisandro Dalcin   PetscInt                np;
69abe9303eSLisandro Dalcin 
70abe9303eSLisandro Dalcin   PetscFunctionBegin;
71abe9303eSLisandro Dalcin   if (p->random) {
72abe9303eSLisandro Dalcin     PetscRandom r;
73abe9303eSLisandro Dalcin     PetscInt   *sizes, *points, v, p;
74abe9303eSLisandro Dalcin     PetscMPIInt rank;
75abe9303eSLisandro Dalcin 
765f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) part), &rank));
775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &r));
785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetInterval(r, 0.0, (PetscScalar) nparts));
795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomSetFromOptions(r));
805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc2(nparts, &sizes, numVertices, &points));
81abe9303eSLisandro Dalcin     for (v = 0; v < numVertices; ++v) {points[v] = v;}
82abe9303eSLisandro Dalcin     for (p = 0; p < nparts; ++p) {sizes[p] = numVertices/nparts + (PetscInt) (p < numVertices % nparts);}
83abe9303eSLisandro Dalcin     for (v = numVertices-1; v > 0; --v) {
84abe9303eSLisandro Dalcin       PetscReal val;
85abe9303eSLisandro Dalcin       PetscInt  w, tmp;
86abe9303eSLisandro Dalcin 
875f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomSetInterval(r, 0.0, (PetscScalar) (v+1)));
885f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscRandomGetValueReal(r, &val));
89abe9303eSLisandro Dalcin       w    = PetscFloorReal(val);
90abe9303eSLisandro Dalcin       tmp       = points[v];
91abe9303eSLisandro Dalcin       points[v] = points[w];
92abe9303eSLisandro Dalcin       points[w] = tmp;
93abe9303eSLisandro Dalcin     }
945f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomDestroy(&r));
955f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPartitionerShellSetPartition(part, nparts, sizes, points));
965f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree2(sizes, points));
97abe9303eSLisandro Dalcin   }
98*28b400f6SJacob Faibussowitsch   PetscCheck(p->section,PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()");
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(p->section, NULL, &np));
1002c71b3e2SJacob Faibussowitsch   PetscCheckFalse(nparts != np,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(ISGetLocalSize(p->partition, &np));
1022c71b3e2SJacob Faibussowitsch   PetscCheckFalse(numVertices != np,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCopy(p->section, partSection));
104abe9303eSLisandro Dalcin   *partition = p->partition;
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectReference((PetscObject) p->partition));
106abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
107abe9303eSLisandro Dalcin }
108abe9303eSLisandro Dalcin 
109abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
110abe9303eSLisandro Dalcin {
111abe9303eSLisandro Dalcin   PetscFunctionBegin;
112abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE; /* PetscPartitionerShell cannot overload the partition call, so it is safe for now */
113abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Shell;
114abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
115abe9303eSLisandro Dalcin   part->ops->reset          = PetscPartitionerReset_Shell;
116abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Shell;
117abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Shell;
118abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
119abe9303eSLisandro Dalcin }
120abe9303eSLisandro Dalcin 
121abe9303eSLisandro Dalcin /*MC
122abe9303eSLisandro Dalcin   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
123abe9303eSLisandro Dalcin 
124abe9303eSLisandro Dalcin   Level: intermediate
125abe9303eSLisandro Dalcin 
126abe9303eSLisandro Dalcin   Options Database Keys:
127abe9303eSLisandro Dalcin .  -petscpartitioner_shell_random - Use a random partition
128abe9303eSLisandro Dalcin 
129abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
130abe9303eSLisandro Dalcin M*/
131abe9303eSLisandro Dalcin 
132abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
133abe9303eSLisandro Dalcin {
134abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p;
135abe9303eSLisandro Dalcin 
136abe9303eSLisandro Dalcin   PetscFunctionBegin;
137abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(part, &p));
139abe9303eSLisandro Dalcin   part->data = p;
140abe9303eSLisandro Dalcin 
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPartitionerInitialize_Shell(part));
142abe9303eSLisandro Dalcin   p->random = PETSC_FALSE;
143abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
144abe9303eSLisandro Dalcin }
145abe9303eSLisandro Dalcin 
146abe9303eSLisandro Dalcin /*@C
147abe9303eSLisandro Dalcin   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
148abe9303eSLisandro Dalcin 
149abe9303eSLisandro Dalcin   Collective on PetscPartitioner
150abe9303eSLisandro Dalcin 
151abe9303eSLisandro Dalcin   Input Parameters:
152abe9303eSLisandro Dalcin + part   - The PetscPartitioner
153abe9303eSLisandro Dalcin . size   - The number of partitions
154abe9303eSLisandro Dalcin . sizes  - array of length size (or NULL) providing the number of points in each partition
155abe9303eSLisandro 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.)
156abe9303eSLisandro Dalcin 
157abe9303eSLisandro Dalcin   Level: developer
158abe9303eSLisandro Dalcin 
159abe9303eSLisandro Dalcin   Notes:
160abe9303eSLisandro Dalcin     It is safe to free the sizes and points arrays after use in this routine.
161abe9303eSLisandro Dalcin 
162abe9303eSLisandro Dalcin .seealso DMPlexDistribute(), PetscPartitionerCreate()
163abe9303eSLisandro Dalcin @*/
164abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[])
165abe9303eSLisandro Dalcin {
166abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
167abe9303eSLisandro Dalcin   PetscInt                proc, numPoints;
168abe9303eSLisandro Dalcin 
169abe9303eSLisandro Dalcin   PetscFunctionBegin;
170abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
171abe9303eSLisandro Dalcin   if (sizes)  {PetscValidPointer(sizes, 3);}
172abe9303eSLisandro Dalcin   if (points) {PetscValidPointer(points, 4);}
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&p->section));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&p->partition));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(p->section, 0, size));
177abe9303eSLisandro Dalcin   if (sizes) {
178abe9303eSLisandro Dalcin     for (proc = 0; proc < size; ++proc) {
1795f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(p->section, proc, sizes[proc]));
180abe9303eSLisandro Dalcin     }
181abe9303eSLisandro Dalcin   }
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(p->section));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(p->section, &numPoints));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition));
185abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
186abe9303eSLisandro Dalcin }
187abe9303eSLisandro Dalcin 
188abe9303eSLisandro Dalcin /*@
189abe9303eSLisandro Dalcin   PetscPartitionerShellSetRandom - Set the flag to use a random partition
190abe9303eSLisandro Dalcin 
191abe9303eSLisandro Dalcin   Collective on PetscPartitioner
192abe9303eSLisandro Dalcin 
193abe9303eSLisandro Dalcin   Input Parameters:
194abe9303eSLisandro Dalcin + part   - The PetscPartitioner
195abe9303eSLisandro Dalcin - random - The flag to use a random partition
196abe9303eSLisandro Dalcin 
197abe9303eSLisandro Dalcin   Level: intermediate
198abe9303eSLisandro Dalcin 
199abe9303eSLisandro Dalcin .seealso PetscPartitionerShellGetRandom(), PetscPartitionerCreate()
200abe9303eSLisandro Dalcin @*/
201abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random)
202abe9303eSLisandro Dalcin {
203abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
204abe9303eSLisandro Dalcin 
205abe9303eSLisandro Dalcin   PetscFunctionBegin;
206abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
207abe9303eSLisandro Dalcin   p->random = random;
208abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
209abe9303eSLisandro Dalcin }
210abe9303eSLisandro Dalcin 
211abe9303eSLisandro Dalcin /*@
212abe9303eSLisandro Dalcin   PetscPartitionerShellGetRandom - get the flag to use a random partition
213abe9303eSLisandro Dalcin 
214abe9303eSLisandro Dalcin   Collective on PetscPartitioner
215abe9303eSLisandro Dalcin 
216abe9303eSLisandro Dalcin   Input Parameter:
217abe9303eSLisandro Dalcin . part   - The PetscPartitioner
218abe9303eSLisandro Dalcin 
219abe9303eSLisandro Dalcin   Output Parameter:
220abe9303eSLisandro Dalcin . random - The flag to use a random partition
221abe9303eSLisandro Dalcin 
222abe9303eSLisandro Dalcin   Level: intermediate
223abe9303eSLisandro Dalcin 
224abe9303eSLisandro Dalcin .seealso PetscPartitionerShellSetRandom(), PetscPartitionerCreate()
225abe9303eSLisandro Dalcin @*/
226abe9303eSLisandro Dalcin PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random)
227abe9303eSLisandro Dalcin {
228abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
229abe9303eSLisandro Dalcin 
230abe9303eSLisandro Dalcin   PetscFunctionBegin;
231abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
232abe9303eSLisandro Dalcin   PetscValidPointer(random, 2);
233abe9303eSLisandro Dalcin   *random = p->random;
234abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
235abe9303eSLisandro Dalcin }
236