xref: /petsc/src/dm/partitioner/impls/shell/partshell.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
99371c9d4SSatish Balay static PetscErrorCode PetscPartitionerReset_Shell(PetscPartitioner part) {
10abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
11abe9303eSLisandro Dalcin 
12abe9303eSLisandro Dalcin   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&p->section));
149566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&p->partition));
15abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
16abe9303eSLisandro Dalcin }
17abe9303eSLisandro Dalcin 
189371c9d4SSatish Balay static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part) {
19abe9303eSLisandro Dalcin   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerReset_Shell(part));
219566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
22abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
23abe9303eSLisandro Dalcin }
24abe9303eSLisandro Dalcin 
259371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Shell_ASCII(PetscPartitioner part, PetscViewer viewer) {
26abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
27abe9303eSLisandro Dalcin 
28abe9303eSLisandro Dalcin   PetscFunctionBegin;
29abe9303eSLisandro Dalcin   if (p->random) {
309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
319566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "using random partition\n"));
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
33abe9303eSLisandro Dalcin   }
34abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
35abe9303eSLisandro Dalcin }
36abe9303eSLisandro Dalcin 
379371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer) {
38abe9303eSLisandro Dalcin   PetscBool iascii;
39abe9303eSLisandro Dalcin 
40abe9303eSLisandro Dalcin   PetscFunctionBegin;
41abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
42abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
439566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
449566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Shell_ASCII(part, viewer));
45abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
46abe9303eSLisandro Dalcin }
47abe9303eSLisandro Dalcin 
489371c9d4SSatish Balay static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscPartitioner part, PetscOptionItems *PetscOptionsObject) {
49d7cc930eSLisandro Dalcin   PetscBool random = PETSC_FALSE, set;
50abe9303eSLisandro Dalcin 
51abe9303eSLisandro Dalcin   PetscFunctionBegin;
52d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner Shell Options");
539566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerShellGetRandom(part, &random));
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &random, &set));
559566063dSJacob Faibussowitsch   if (set) PetscCall(PetscPartitionerShellSetRandom(part, random));
56d0609cedSBarry Smith   PetscOptionsHeadEnd();
57abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
58abe9303eSLisandro Dalcin }
59abe9303eSLisandro Dalcin 
609371c9d4SSatish Balay static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) {
61abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
62abe9303eSLisandro Dalcin   PetscInt                np;
63abe9303eSLisandro Dalcin 
64abe9303eSLisandro Dalcin   PetscFunctionBegin;
65abe9303eSLisandro Dalcin   if (p->random) {
66abe9303eSLisandro Dalcin     PetscRandom r;
67abe9303eSLisandro Dalcin     PetscInt   *sizes, *points, v, p;
68abe9303eSLisandro Dalcin     PetscMPIInt rank;
69abe9303eSLisandro Dalcin 
709566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank));
719566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
729566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)nparts));
739566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(r));
749566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nparts, &sizes, numVertices, &points));
75abe9303eSLisandro Dalcin     for (v = 0; v < numVertices; ++v) { points[v] = v; }
76abe9303eSLisandro Dalcin     for (p = 0; p < nparts; ++p) { sizes[p] = numVertices / nparts + (PetscInt)(p < numVertices % nparts); }
77abe9303eSLisandro Dalcin     for (v = numVertices - 1; v > 0; --v) {
78abe9303eSLisandro Dalcin       PetscReal val;
79abe9303eSLisandro Dalcin       PetscInt  w, tmp;
80abe9303eSLisandro Dalcin 
819566063dSJacob Faibussowitsch       PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)(v + 1)));
829566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(r, &val));
83abe9303eSLisandro Dalcin       w         = PetscFloorReal(val);
84abe9303eSLisandro Dalcin       tmp       = points[v];
85abe9303eSLisandro Dalcin       points[v] = points[w];
86abe9303eSLisandro Dalcin       points[w] = tmp;
87abe9303eSLisandro Dalcin     }
889566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&r));
899566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, nparts, sizes, points));
909566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
91abe9303eSLisandro Dalcin   }
9228b400f6SJacob Faibussowitsch   PetscCheck(p->section, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()");
939566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(p->section, NULL, &np));
9463a3b9bcSJacob Faibussowitsch   PetscCheck(nparts == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %" PetscInt_FMT " != configured partitions %" PetscInt_FMT, nparts, np);
959566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(p->partition, &np));
9663a3b9bcSJacob Faibussowitsch   PetscCheck(numVertices == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %" PetscInt_FMT " != configured vertices %" PetscInt_FMT, numVertices, np);
979566063dSJacob Faibussowitsch   PetscCall(PetscSectionCopy(p->section, partSection));
98abe9303eSLisandro Dalcin   *partition = p->partition;
999566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)p->partition));
100abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
101abe9303eSLisandro Dalcin }
102abe9303eSLisandro Dalcin 
1039371c9d4SSatish Balay static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part) {
104abe9303eSLisandro Dalcin   PetscFunctionBegin;
105abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE; /* PetscPartitionerShell cannot overload the partition call, so it is safe for now */
106abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Shell;
107abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
108abe9303eSLisandro Dalcin   part->ops->reset          = PetscPartitionerReset_Shell;
109abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Shell;
110abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Shell;
111abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
112abe9303eSLisandro Dalcin }
113abe9303eSLisandro Dalcin 
114abe9303eSLisandro Dalcin /*MC
115abe9303eSLisandro Dalcin   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
116abe9303eSLisandro Dalcin 
117abe9303eSLisandro Dalcin   Level: intermediate
118abe9303eSLisandro Dalcin 
119abe9303eSLisandro Dalcin   Options Database Keys:
120abe9303eSLisandro Dalcin .  -petscpartitioner_shell_random - Use a random partition
121abe9303eSLisandro Dalcin 
122db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
123abe9303eSLisandro Dalcin M*/
124abe9303eSLisandro Dalcin 
1259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part) {
126abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p;
127abe9303eSLisandro Dalcin 
128abe9303eSLisandro Dalcin   PetscFunctionBegin;
129abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1309566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(part, &p));
131abe9303eSLisandro Dalcin   part->data = p;
132abe9303eSLisandro Dalcin 
1339566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Shell(part));
134abe9303eSLisandro Dalcin   p->random = PETSC_FALSE;
135abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
136abe9303eSLisandro Dalcin }
137abe9303eSLisandro Dalcin 
138abe9303eSLisandro Dalcin /*@C
139abe9303eSLisandro Dalcin   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
140abe9303eSLisandro Dalcin 
141abe9303eSLisandro Dalcin   Collective on PetscPartitioner
142abe9303eSLisandro Dalcin 
143abe9303eSLisandro Dalcin   Input Parameters:
144abe9303eSLisandro Dalcin + part   - The PetscPartitioner
145abe9303eSLisandro Dalcin . size   - The number of partitions
146abe9303eSLisandro Dalcin . sizes  - array of length size (or NULL) providing the number of points in each partition
147abe9303eSLisandro 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.)
148abe9303eSLisandro Dalcin 
149abe9303eSLisandro Dalcin   Level: developer
150abe9303eSLisandro Dalcin 
151abe9303eSLisandro Dalcin   Notes:
152abe9303eSLisandro Dalcin     It is safe to free the sizes and points arrays after use in this routine.
153abe9303eSLisandro Dalcin 
154db781477SPatrick Sanan .seealso `DMPlexDistribute()`, `PetscPartitionerCreate()`
155abe9303eSLisandro Dalcin @*/
1569371c9d4SSatish Balay PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[]) {
157abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
158abe9303eSLisandro Dalcin   PetscInt                proc, numPoints;
159abe9303eSLisandro Dalcin 
160abe9303eSLisandro Dalcin   PetscFunctionBegin;
161abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
162dadcf809SJacob Faibussowitsch   if (sizes) { PetscValidIntPointer(sizes, 3); }
163dadcf809SJacob Faibussowitsch   if (points) { PetscValidIntPointer(points, 4); }
1649566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&p->section));
1659566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&p->partition));
1669566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)part), &p->section));
1679566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(p->section, 0, size));
168abe9303eSLisandro Dalcin   if (sizes) {
169*48a46eb9SPierre Jolivet     for (proc = 0; proc < size; ++proc) PetscCall(PetscSectionSetDof(p->section, proc, sizes[proc]));
170abe9303eSLisandro Dalcin   }
1719566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(p->section));
1729566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(p->section, &numPoints));
1739566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), numPoints, points, PETSC_COPY_VALUES, &p->partition));
174abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
175abe9303eSLisandro Dalcin }
176abe9303eSLisandro Dalcin 
177abe9303eSLisandro Dalcin /*@
178abe9303eSLisandro Dalcin   PetscPartitionerShellSetRandom - Set the flag to use a random partition
179abe9303eSLisandro Dalcin 
180abe9303eSLisandro Dalcin   Collective on PetscPartitioner
181abe9303eSLisandro Dalcin 
182abe9303eSLisandro Dalcin   Input Parameters:
183abe9303eSLisandro Dalcin + part   - The PetscPartitioner
184abe9303eSLisandro Dalcin - random - The flag to use a random partition
185abe9303eSLisandro Dalcin 
186abe9303eSLisandro Dalcin   Level: intermediate
187abe9303eSLisandro Dalcin 
188db781477SPatrick Sanan .seealso `PetscPartitionerShellGetRandom()`, `PetscPartitionerCreate()`
189abe9303eSLisandro Dalcin @*/
1909371c9d4SSatish Balay PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random) {
191abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
192abe9303eSLisandro Dalcin 
193abe9303eSLisandro Dalcin   PetscFunctionBegin;
194abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
195abe9303eSLisandro Dalcin   p->random = random;
196abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
197abe9303eSLisandro Dalcin }
198abe9303eSLisandro Dalcin 
199abe9303eSLisandro Dalcin /*@
200abe9303eSLisandro Dalcin   PetscPartitionerShellGetRandom - get the flag to use a random partition
201abe9303eSLisandro Dalcin 
202abe9303eSLisandro Dalcin   Collective on PetscPartitioner
203abe9303eSLisandro Dalcin 
204abe9303eSLisandro Dalcin   Input Parameter:
205abe9303eSLisandro Dalcin . part   - The PetscPartitioner
206abe9303eSLisandro Dalcin 
207abe9303eSLisandro Dalcin   Output Parameter:
208abe9303eSLisandro Dalcin . random - The flag to use a random partition
209abe9303eSLisandro Dalcin 
210abe9303eSLisandro Dalcin   Level: intermediate
211abe9303eSLisandro Dalcin 
212db781477SPatrick Sanan .seealso `PetscPartitionerShellSetRandom()`, `PetscPartitionerCreate()`
213abe9303eSLisandro Dalcin @*/
2149371c9d4SSatish Balay PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random) {
215abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
216abe9303eSLisandro Dalcin 
217abe9303eSLisandro Dalcin   PetscFunctionBegin;
218abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
219dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(random, 2);
220abe9303eSLisandro Dalcin   *random = p->random;
221abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
222abe9303eSLisandro Dalcin }
223