xref: /petsc/src/dm/partitioner/impls/shell/partshell.c (revision bcda9346efad4e5ba2d553af84eb238771ba1e25)
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 
PetscPartitionerReset_Shell(PetscPartitioner part)9d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerReset_Shell(PetscPartitioner part)
10d71ae5a4SJacob Faibussowitsch {
11abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
12abe9303eSLisandro Dalcin 
13abe9303eSLisandro Dalcin   PetscFunctionBegin;
149566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&p->section));
159566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&p->partition));
163ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17abe9303eSLisandro Dalcin }
18abe9303eSLisandro Dalcin 
PetscPartitionerDestroy_Shell(PetscPartitioner part)19d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
20d71ae5a4SJacob Faibussowitsch {
21abe9303eSLisandro Dalcin   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerReset_Shell(part));
239566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25abe9303eSLisandro Dalcin }
26abe9303eSLisandro Dalcin 
PetscPartitionerView_Shell_ASCII(PetscPartitioner part,PetscViewer viewer)27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Shell_ASCII(PetscPartitioner part, PetscViewer viewer)
28d71ae5a4SJacob Faibussowitsch {
29abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
30abe9303eSLisandro Dalcin 
31abe9303eSLisandro Dalcin   PetscFunctionBegin;
32abe9303eSLisandro Dalcin   if (p->random) {
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPushTab(viewer));
349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "using random partition\n"));
359566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
36abe9303eSLisandro Dalcin   }
373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38abe9303eSLisandro Dalcin }
39abe9303eSLisandro Dalcin 
PetscPartitionerView_Shell(PetscPartitioner part,PetscViewer viewer)40d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
41d71ae5a4SJacob Faibussowitsch {
42*9f196a02SMartin Diehl   PetscBool isascii;
43abe9303eSLisandro Dalcin 
44abe9303eSLisandro Dalcin   PetscFunctionBegin;
45abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
46abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
47*9f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
48*9f196a02SMartin Diehl   if (isascii) PetscCall(PetscPartitionerView_Shell_ASCII(part, viewer));
493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
50abe9303eSLisandro Dalcin }
51abe9303eSLisandro Dalcin 
PetscPartitionerSetFromOptions_Shell(PetscPartitioner part,PetscOptionItems PetscOptionsObject)52ce78bad3SBarry Smith static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscPartitioner part, PetscOptionItems PetscOptionsObject)
53d71ae5a4SJacob Faibussowitsch {
546166d6f7SMatthew G. Knepley   PetscInt    sizes[16], points[1024];
556166d6f7SMatthew G. Knepley   PetscInt    Npart = 16, Npoints = 1024;
566166d6f7SMatthew G. Knepley   PetscBool   random = PETSC_FALSE, set, flgSizes, flgPoints;
576166d6f7SMatthew G. Knepley   PetscMPIInt rank;
58abe9303eSLisandro Dalcin 
59abe9303eSLisandro Dalcin   PetscFunctionBegin;
606166d6f7SMatthew G. Knepley   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank));
61d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner Shell Options");
629566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerShellGetRandom(part, &random));
639566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &random, &set));
649566063dSJacob Faibussowitsch   if (set) PetscCall(PetscPartitionerShellSetRandom(part, random));
656166d6f7SMatthew G. Knepley   PetscCall(PetscOptionsIntArray("-petscpartitioner_shell_sizes", "The size of each partition on rank 0", "PetscPartitionerShellSetPartition", sizes, &Npart, &flgSizes));
666166d6f7SMatthew G. Knepley   PetscCall(PetscOptionsIntArray("-petscpartitioner_shell_points", "The points in each partition on rank 0", "PetscPartitionerShellSetPartition", points, &Npoints, &flgPoints));
676166d6f7SMatthew G. Knepley   PetscCheck(!(flgSizes ^ flgPoints), PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_WRONG, "Must specify both the partition sizes and points");
686166d6f7SMatthew G. Knepley   if (flgSizes) {
696166d6f7SMatthew G. Knepley     PetscInt Np = 0;
706166d6f7SMatthew G. Knepley 
716166d6f7SMatthew G. Knepley     for (PetscInt i = 0; i < Npart; ++i) Np += sizes[i];
726166d6f7SMatthew G. Knepley     PetscCheck(Np == Npoints, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_WRONG, "Number of input points %" PetscInt_FMT " != %" PetscInt_FMT " sum of partition sizes", Npoints, Np);
736166d6f7SMatthew G. Knepley     if (!rank) PetscCall(PetscPartitionerShellSetPartition(part, Npart, sizes, points));
746166d6f7SMatthew G. Knepley     else {
756166d6f7SMatthew G. Knepley       PetscCall(PetscArrayzero(sizes, Npart));
766166d6f7SMatthew G. Knepley       PetscCall(PetscPartitionerShellSetPartition(part, Npart, sizes, points));
776166d6f7SMatthew G. Knepley     }
786166d6f7SMatthew G. Knepley   }
79d0609cedSBarry Smith   PetscOptionsHeadEnd();
803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
81abe9303eSLisandro Dalcin }
82abe9303eSLisandro Dalcin 
PetscPartitionerPartition_Shell(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * partition)8321c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
84d71ae5a4SJacob Faibussowitsch {
85abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
86abe9303eSLisandro Dalcin   PetscInt                np;
87abe9303eSLisandro Dalcin 
88abe9303eSLisandro Dalcin   PetscFunctionBegin;
89abe9303eSLisandro Dalcin   if (p->random) {
90abe9303eSLisandro Dalcin     PetscRandom r;
91abe9303eSLisandro Dalcin     PetscInt   *sizes, *points, v, p;
92abe9303eSLisandro Dalcin     PetscMPIInt rank;
93abe9303eSLisandro Dalcin 
949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank));
959566063dSJacob Faibussowitsch     PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
969566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)nparts));
979566063dSJacob Faibussowitsch     PetscCall(PetscRandomSetFromOptions(r));
989566063dSJacob Faibussowitsch     PetscCall(PetscCalloc2(nparts, &sizes, numVertices, &points));
99ad540459SPierre Jolivet     for (v = 0; v < numVertices; ++v) points[v] = v;
100ad540459SPierre Jolivet     for (p = 0; p < nparts; ++p) sizes[p] = numVertices / nparts + (PetscInt)(p < numVertices % nparts);
101abe9303eSLisandro Dalcin     for (v = numVertices - 1; v > 0; --v) {
102abe9303eSLisandro Dalcin       PetscReal val;
103abe9303eSLisandro Dalcin       PetscInt  w, tmp;
104abe9303eSLisandro Dalcin 
1059566063dSJacob Faibussowitsch       PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)(v + 1)));
1069566063dSJacob Faibussowitsch       PetscCall(PetscRandomGetValueReal(r, &val));
107abe9303eSLisandro Dalcin       w         = PetscFloorReal(val);
108abe9303eSLisandro Dalcin       tmp       = points[v];
109abe9303eSLisandro Dalcin       points[v] = points[w];
110abe9303eSLisandro Dalcin       points[w] = tmp;
111abe9303eSLisandro Dalcin     }
1129566063dSJacob Faibussowitsch     PetscCall(PetscRandomDestroy(&r));
1139566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, nparts, sizes, points));
1149566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
115abe9303eSLisandro Dalcin   }
11628b400f6SJacob Faibussowitsch   PetscCheck(p->section, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()");
1179566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(p->section, NULL, &np));
11863a3b9bcSJacob Faibussowitsch   PetscCheck(nparts == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %" PetscInt_FMT " != configured partitions %" PetscInt_FMT, nparts, np);
1199566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(p->partition, &np));
12063a3b9bcSJacob Faibussowitsch   PetscCheck(numVertices == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %" PetscInt_FMT " != configured vertices %" PetscInt_FMT, numVertices, np);
1219566063dSJacob Faibussowitsch   PetscCall(PetscSectionCopy(p->section, partSection));
122abe9303eSLisandro Dalcin   *partition = p->partition;
1239566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)p->partition));
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
125abe9303eSLisandro Dalcin }
126abe9303eSLisandro Dalcin 
PetscPartitionerInitialize_Shell(PetscPartitioner part)127d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
128d71ae5a4SJacob Faibussowitsch {
129abe9303eSLisandro Dalcin   PetscFunctionBegin;
130abe9303eSLisandro Dalcin   part->noGraph             = PETSC_TRUE; /* PetscPartitionerShell cannot overload the partition call, so it is safe for now */
131abe9303eSLisandro Dalcin   part->ops->view           = PetscPartitionerView_Shell;
132abe9303eSLisandro Dalcin   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
133abe9303eSLisandro Dalcin   part->ops->reset          = PetscPartitionerReset_Shell;
134abe9303eSLisandro Dalcin   part->ops->destroy        = PetscPartitionerDestroy_Shell;
135abe9303eSLisandro Dalcin   part->ops->partition      = PetscPartitionerPartition_Shell;
1363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
137abe9303eSLisandro Dalcin }
138abe9303eSLisandro Dalcin 
139abe9303eSLisandro Dalcin /*MC
140abe9303eSLisandro Dalcin   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
141abe9303eSLisandro Dalcin 
142abe9303eSLisandro Dalcin   Level: intermediate
143abe9303eSLisandro Dalcin 
144abe9303eSLisandro Dalcin   Options Database Keys:
145abe9303eSLisandro Dalcin .  -petscpartitioner_shell_random - Use a random partition
146abe9303eSLisandro Dalcin 
147db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
148abe9303eSLisandro Dalcin M*/
149abe9303eSLisandro Dalcin 
PetscPartitionerCreate_Shell(PetscPartitioner part)150d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
151d71ae5a4SJacob Faibussowitsch {
152abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p;
153abe9303eSLisandro Dalcin 
154abe9303eSLisandro Dalcin   PetscFunctionBegin;
155abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1564dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
157abe9303eSLisandro Dalcin   part->data = p;
158abe9303eSLisandro Dalcin 
1599566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Shell(part));
160abe9303eSLisandro Dalcin   p->random = PETSC_FALSE;
1613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162abe9303eSLisandro Dalcin }
163abe9303eSLisandro Dalcin 
164cc4c1da9SBarry Smith /*@
16535cb6cd3SPierre Jolivet   PetscPartitionerShellSetPartition - Set an artificial partition for a mesh
166abe9303eSLisandro Dalcin 
16720f4b53cSBarry Smith   Collective
168abe9303eSLisandro Dalcin 
169abe9303eSLisandro Dalcin   Input Parameters:
17020f4b53cSBarry Smith + part   - The `PetscPartitioner`
171abe9303eSLisandro Dalcin . size   - The number of partitions
17220f4b53cSBarry Smith . sizes  - array of length size (or `NULL`) providing the number of points in each partition
17320f4b53cSBarry Smith - 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.)
174abe9303eSLisandro Dalcin 
175abe9303eSLisandro Dalcin   Level: developer
176abe9303eSLisandro Dalcin 
17720f4b53cSBarry Smith   Note:
178abe9303eSLisandro Dalcin   It is safe to free the sizes and points arrays after use in this routine.
179abe9303eSLisandro Dalcin 
180a4e35b19SJacob Faibussowitsch .seealso: `DMPlexDistribute()`, `PetscPartitionerCreate()`
181abe9303eSLisandro Dalcin @*/
PetscPartitionerShellSetPartition(PetscPartitioner part,PetscInt size,const PetscInt sizes[],const PetscInt points[])182d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[])
183d71ae5a4SJacob Faibussowitsch {
184abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
185abe9303eSLisandro Dalcin   PetscInt                proc, numPoints;
186abe9303eSLisandro Dalcin 
187abe9303eSLisandro Dalcin   PetscFunctionBegin;
188abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
1894f572ea9SToby Isaac   if (sizes) PetscAssertPointer(sizes, 3);
1904f572ea9SToby Isaac   if (points) PetscAssertPointer(points, 4);
1919566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&p->section));
1929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&p->partition));
1939566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)part), &p->section));
1949566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(p->section, 0, size));
195abe9303eSLisandro Dalcin   if (sizes) {
19648a46eb9SPierre Jolivet     for (proc = 0; proc < size; ++proc) PetscCall(PetscSectionSetDof(p->section, proc, sizes[proc]));
197abe9303eSLisandro Dalcin   }
1989566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(p->section));
1999566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(p->section, &numPoints));
2009566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), numPoints, points, PETSC_COPY_VALUES, &p->partition));
2013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
202abe9303eSLisandro Dalcin }
203abe9303eSLisandro Dalcin 
204abe9303eSLisandro Dalcin /*@
205abe9303eSLisandro Dalcin   PetscPartitionerShellSetRandom - Set the flag to use a random partition
206abe9303eSLisandro Dalcin 
20720f4b53cSBarry Smith   Collective
208abe9303eSLisandro Dalcin 
209abe9303eSLisandro Dalcin   Input Parameters:
21020f4b53cSBarry Smith + part   - The `PetscPartitioner`
211abe9303eSLisandro Dalcin - random - The flag to use a random partition
212abe9303eSLisandro Dalcin 
213abe9303eSLisandro Dalcin   Level: intermediate
214abe9303eSLisandro Dalcin 
215a4e35b19SJacob Faibussowitsch .seealso: `PetscPartitionerShellGetRandom()`, `PetscPartitionerCreate()`
216abe9303eSLisandro Dalcin @*/
PetscPartitionerShellSetRandom(PetscPartitioner part,PetscBool random)217d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random)
218d71ae5a4SJacob Faibussowitsch {
219abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
220abe9303eSLisandro Dalcin 
221abe9303eSLisandro Dalcin   PetscFunctionBegin;
222abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
223abe9303eSLisandro Dalcin   p->random = random;
2243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
225abe9303eSLisandro Dalcin }
226abe9303eSLisandro Dalcin 
227abe9303eSLisandro Dalcin /*@
228abe9303eSLisandro Dalcin   PetscPartitionerShellGetRandom - get the flag to use a random partition
229abe9303eSLisandro Dalcin 
23020f4b53cSBarry Smith   Collective
231abe9303eSLisandro Dalcin 
232abe9303eSLisandro Dalcin   Input Parameter:
23320f4b53cSBarry Smith . part - The `PetscPartitioner`
234abe9303eSLisandro Dalcin 
235abe9303eSLisandro Dalcin   Output Parameter:
236abe9303eSLisandro Dalcin . random - The flag to use a random partition
237abe9303eSLisandro Dalcin 
238abe9303eSLisandro Dalcin   Level: intermediate
239abe9303eSLisandro Dalcin 
240a4e35b19SJacob Faibussowitsch .seealso: `PetscPartitionerShellSetRandom()`, `PetscPartitionerCreate()`
241abe9303eSLisandro Dalcin @*/
PetscPartitionerShellGetRandom(PetscPartitioner part,PetscBool * random)242d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random)
243d71ae5a4SJacob Faibussowitsch {
244abe9303eSLisandro Dalcin   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data;
245abe9303eSLisandro Dalcin 
246abe9303eSLisandro Dalcin   PetscFunctionBegin;
247abe9303eSLisandro Dalcin   PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL);
2484f572ea9SToby Isaac   PetscAssertPointer(random, 2);
249abe9303eSLisandro Dalcin   *random = p->random;
2503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
251abe9303eSLisandro Dalcin }
252