1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
2abe9303eSLisandro Dalcin
3cc4c1da9SBarry Smith /*@
420f4b53cSBarry Smith PetscPartitionerSetType - Builds a particular `PetscPartitioner`
5abe9303eSLisandro Dalcin
620f4b53cSBarry Smith Collective
7abe9303eSLisandro Dalcin
8abe9303eSLisandro Dalcin Input Parameters:
920f4b53cSBarry Smith + part - The `PetscPartitioner` object
10abe9303eSLisandro Dalcin - name - The kind of partitioner
11abe9303eSLisandro Dalcin
12abe9303eSLisandro Dalcin Options Database Key:
1320f4b53cSBarry Smith . -petscpartitioner_type <type> - Sets the `PetscPartitioner` type
146e0dc72eSMatthew G. Knepley
15abe9303eSLisandro Dalcin Level: intermediate
16abe9303eSLisandro Dalcin
1720f4b53cSBarry Smith Note:
1820f4b53cSBarry Smith .vb
1920f4b53cSBarry Smith PETSCPARTITIONERCHACO - The Chaco partitioner (--download-chaco)
2020f4b53cSBarry Smith PETSCPARTITIONERPARMETIS - The ParMetis partitioner (--download-parmetis)
2120f4b53cSBarry Smith PETSCPARTITIONERSHELL - A shell partitioner implemented by the user
2220f4b53cSBarry Smith PETSCPARTITIONERSIMPLE - A simple partitioner that divides cells into equal, contiguous chunks
2320f4b53cSBarry Smith PETSCPARTITIONERGATHER - Gathers all cells onto process 0
2420f4b53cSBarry Smith .ve
2520f4b53cSBarry Smith
26db781477SPatrick Sanan .seealso: `PetscPartitionerGetType()`, `PetscPartitionerCreate()`
27abe9303eSLisandro Dalcin @*/
PetscPartitionerSetType(PetscPartitioner part,PetscPartitionerType name)28d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
29d71ae5a4SJacob Faibussowitsch {
30abe9303eSLisandro Dalcin PetscErrorCode (*r)(PetscPartitioner);
31abe9303eSLisandro Dalcin PetscBool match;
32abe9303eSLisandro Dalcin
33abe9303eSLisandro Dalcin PetscFunctionBegin;
34abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)part, name, &match));
363ba16761SJacob Faibussowitsch if (match) PetscFunctionReturn(PETSC_SUCCESS);
37abe9303eSLisandro Dalcin
389566063dSJacob Faibussowitsch PetscCall(PetscPartitionerRegisterAll());
399566063dSJacob Faibussowitsch PetscCall(PetscFunctionListFind(PetscPartitionerList, name, &r));
4028b400f6SJacob Faibussowitsch PetscCheck(r, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);
41abe9303eSLisandro Dalcin
42dbbe0bcdSBarry Smith PetscTryTypeMethod(part, destroy);
43abe9303eSLisandro Dalcin part->noGraph = PETSC_FALSE;
449566063dSJacob Faibussowitsch PetscCall(PetscMemzero(part->ops, sizeof(*part->ops)));
459566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)part, name));
469566063dSJacob Faibussowitsch PetscCall((*r)(part));
473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
48abe9303eSLisandro Dalcin }
49abe9303eSLisandro Dalcin
50cc4c1da9SBarry Smith /*@
51abe9303eSLisandro Dalcin PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.
52abe9303eSLisandro Dalcin
53abe9303eSLisandro Dalcin Not Collective
54abe9303eSLisandro Dalcin
55abe9303eSLisandro Dalcin Input Parameter:
56abe9303eSLisandro Dalcin . part - The PetscPartitioner
57abe9303eSLisandro Dalcin
58abe9303eSLisandro Dalcin Output Parameter:
59abe9303eSLisandro Dalcin . name - The PetscPartitioner type name
60abe9303eSLisandro Dalcin
61abe9303eSLisandro Dalcin Level: intermediate
62abe9303eSLisandro Dalcin
63db781477SPatrick Sanan .seealso: `PetscPartitionerSetType()`, `PetscPartitionerCreate()`
64abe9303eSLisandro Dalcin @*/
PetscPartitionerGetType(PetscPartitioner part,PetscPartitionerType * name)65d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
66d71ae5a4SJacob Faibussowitsch {
67abe9303eSLisandro Dalcin PetscFunctionBegin;
68abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
694f572ea9SToby Isaac PetscAssertPointer(name, 2);
70abe9303eSLisandro Dalcin *name = ((PetscObject)part)->type_name;
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72abe9303eSLisandro Dalcin }
73abe9303eSLisandro Dalcin
74ffeef943SBarry Smith /*@
7520f4b53cSBarry Smith PetscPartitionerViewFromOptions - View a `PetscPartitioner` object based on options in the options database
76abe9303eSLisandro Dalcin
7720f4b53cSBarry Smith Collective
78abe9303eSLisandro Dalcin
79abe9303eSLisandro Dalcin Input Parameters:
8020f4b53cSBarry Smith + A - the `PetscPartitioner` object
8120f4b53cSBarry Smith . obj - Optional `PetscObject` that provides the options prefix
82abe9303eSLisandro Dalcin - name - command line option
83abe9303eSLisandro Dalcin
84abe9303eSLisandro Dalcin Level: intermediate
8520f4b53cSBarry Smith
8620f4b53cSBarry Smith Note:
8720f4b53cSBarry Smith See `PetscObjectViewFromOptions()` for the various forms of viewers that may be used
8820f4b53cSBarry Smith
89db781477SPatrick Sanan .seealso: `PetscPartitionerView()`, `PetscObjectViewFromOptions()`
90abe9303eSLisandro Dalcin @*/
PetscPartitionerViewFromOptions(PetscPartitioner A,PetscObject obj,const char name[])91d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A, PetscObject obj, const char name[])
92d71ae5a4SJacob Faibussowitsch {
93abe9303eSLisandro Dalcin PetscFunctionBegin;
94abe9303eSLisandro Dalcin PetscValidHeaderSpecific(A, PETSCPARTITIONER_CLASSID, 1);
959566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97abe9303eSLisandro Dalcin }
98abe9303eSLisandro Dalcin
99abe9303eSLisandro Dalcin /*@
10020f4b53cSBarry Smith PetscPartitionerView - Views a `PetscPartitioner`
101abe9303eSLisandro Dalcin
10220f4b53cSBarry Smith Collective
103abe9303eSLisandro Dalcin
104d8d19677SJose E. Roman Input Parameters:
10520f4b53cSBarry Smith + part - the `PetscPartitioner` object to view
106abe9303eSLisandro Dalcin - v - the viewer
107abe9303eSLisandro Dalcin
108abe9303eSLisandro Dalcin Level: developer
109abe9303eSLisandro Dalcin
110db781477SPatrick Sanan .seealso: `PetscPartitionerDestroy()`
111abe9303eSLisandro Dalcin @*/
PetscPartitionerView(PetscPartitioner part,PetscViewer v)112d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
113d71ae5a4SJacob Faibussowitsch {
114abe9303eSLisandro Dalcin PetscMPIInt size;
115abe9303eSLisandro Dalcin PetscBool isascii;
116abe9303eSLisandro Dalcin
117abe9303eSLisandro Dalcin PetscFunctionBegin;
118abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1199566063dSJacob Faibussowitsch if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)part), &v));
1209566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
121*3b9d9b65SStefano Zampini if (isascii && part->printHeader) {
1229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)part), &size));
1239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, "Graph Partitioner: %d MPI Process%s\n", size, size > 1 ? "es" : ""));
1249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " type: %s\n", ((PetscObject)part)->type_name));
12563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " edge cut: %" PetscInt_FMT "\n", part->edgeCut));
12663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " balance: %.2g\n", (double)part->balance));
1279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(v, " use vertex weights: %d\n", part->usevwgt));
12821c92275SMatthew G. Knepley PetscCall(PetscViewerASCIIPrintf(v, " use edge weights: %d\n", part->useewgt));
129abe9303eSLisandro Dalcin }
130dbbe0bcdSBarry Smith PetscTryTypeMethod(part, view, v);
1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
132abe9303eSLisandro Dalcin }
133abe9303eSLisandro Dalcin
PetscPartitionerGetDefaultType(MPI_Comm comm,const char ** defaultType)134d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerGetDefaultType(MPI_Comm comm, const char **defaultType)
135d71ae5a4SJacob Faibussowitsch {
136abe9303eSLisandro Dalcin PetscMPIInt size;
137abe9303eSLisandro Dalcin
138abe9303eSLisandro Dalcin PetscFunctionBegin;
1399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
140abe9303eSLisandro Dalcin if (size == 1) {
141abe9303eSLisandro Dalcin *defaultType = PETSCPARTITIONERSIMPLE;
142abe9303eSLisandro Dalcin } else {
143abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_PARMETIS)
144abe9303eSLisandro Dalcin *defaultType = PETSCPARTITIONERPARMETIS;
145abe9303eSLisandro Dalcin #elif defined(PETSC_HAVE_PTSCOTCH)
146abe9303eSLisandro Dalcin *defaultType = PETSCPARTITIONERPTSCOTCH;
147abe9303eSLisandro Dalcin #elif defined(PETSC_HAVE_CHACO)
148abe9303eSLisandro Dalcin *defaultType = PETSCPARTITIONERCHACO;
149abe9303eSLisandro Dalcin #else
150abe9303eSLisandro Dalcin *defaultType = PETSCPARTITIONERSIMPLE;
151abe9303eSLisandro Dalcin #endif
152abe9303eSLisandro Dalcin }
1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
154abe9303eSLisandro Dalcin }
155abe9303eSLisandro Dalcin
156abe9303eSLisandro Dalcin /*@
15720f4b53cSBarry Smith PetscPartitionerSetFromOptions - sets parameters in a `PetscPartitioner` from the options database
158abe9303eSLisandro Dalcin
15920f4b53cSBarry Smith Collective
160abe9303eSLisandro Dalcin
161abe9303eSLisandro Dalcin Input Parameter:
16220f4b53cSBarry Smith . part - the `PetscPartitioner` object to set options for
163abe9303eSLisandro Dalcin
164abe9303eSLisandro Dalcin Options Database Keys:
16520f4b53cSBarry Smith + -petscpartitioner_type <type> - Sets the `PetscPartitioner` type; use -help for a list of available types
166abe9303eSLisandro Dalcin . -petscpartitioner_use_vertex_weights - Uses weights associated with the graph vertices
167648c30bcSBarry Smith - -petscpartitioner_view_graph - View the graph each time PetscPartitionerPartition is called. Viewer can be customized, see `PetscOptionsCreateViewer()`
168abe9303eSLisandro Dalcin
169abe9303eSLisandro Dalcin Level: developer
170abe9303eSLisandro Dalcin
171db781477SPatrick Sanan .seealso: `PetscPartitionerView()`, `PetscPartitionerSetType()`, `PetscPartitionerPartition()`
172abe9303eSLisandro Dalcin @*/
PetscPartitionerSetFromOptions(PetscPartitioner part)173d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
174d71ae5a4SJacob Faibussowitsch {
175d7cc930eSLisandro Dalcin const char *currentType = NULL;
176abe9303eSLisandro Dalcin char name[256];
177abe9303eSLisandro Dalcin PetscBool flg;
178abe9303eSLisandro Dalcin
179abe9303eSLisandro Dalcin PetscFunctionBegin;
180abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
181d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)part);
1829566063dSJacob Faibussowitsch PetscCall(PetscPartitionerGetType(part, ¤tType));
1839566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, currentType, name, sizeof(name), &flg));
1841baa6e33SBarry Smith if (flg) PetscCall(PetscPartitionerSetType(part, name));
1859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscpartitioner_use_vertex_weights", "Use vertex weights", "", part->usevwgt, &part->usevwgt, NULL));
18621c92275SMatthew G. Knepley PetscCall(PetscOptionsBool("-petscpartitioner_use_edge_weights", "Use edge weights", "", part->useewgt, &part->useewgt, NULL));
187dbbe0bcdSBarry Smith PetscTryTypeMethod(part, setfromoptions, PetscOptionsObject);
188648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&part->viewer));
189648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&part->viewerGraph));
190*3b9d9b65SStefano Zampini PetscCall(PetscOptionsCreateViewer(((PetscObject)part)->comm, ((PetscObject)part)->options, ((PetscObject)part)->prefix, "-petscpartitioner_view", &part->viewer, &part->viewerFmt, NULL));
191*3b9d9b65SStefano Zampini PetscCall(PetscOptionsCreateViewer(((PetscObject)part)->comm, ((PetscObject)part)->options, ((PetscObject)part)->prefix, "-petscpartitioner_view_graph", &part->viewerGraph, NULL, NULL));
192abe9303eSLisandro Dalcin /* process any options handlers added with PetscObjectAddOptionsHandler() */
193dbbe0bcdSBarry Smith PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)part, PetscOptionsObject));
194d0609cedSBarry Smith PetscOptionsEnd();
1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
196abe9303eSLisandro Dalcin }
197abe9303eSLisandro Dalcin
198abe9303eSLisandro Dalcin /*@
19920f4b53cSBarry Smith PetscPartitionerSetUp - Construct data structures for the `PetscPartitioner`
200abe9303eSLisandro Dalcin
20120f4b53cSBarry Smith Collective
202abe9303eSLisandro Dalcin
203abe9303eSLisandro Dalcin Input Parameter:
20420f4b53cSBarry Smith . part - the `PetscPartitioner` object to setup
205abe9303eSLisandro Dalcin
206abe9303eSLisandro Dalcin Level: developer
207abe9303eSLisandro Dalcin
208db781477SPatrick Sanan .seealso: `PetscPartitionerView()`, `PetscPartitionerDestroy()`
209abe9303eSLisandro Dalcin @*/
PetscPartitionerSetUp(PetscPartitioner part)210d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
211d71ae5a4SJacob Faibussowitsch {
212abe9303eSLisandro Dalcin PetscFunctionBegin;
213abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
214dbbe0bcdSBarry Smith PetscTryTypeMethod(part, setup);
2153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
216abe9303eSLisandro Dalcin }
217abe9303eSLisandro Dalcin
218abe9303eSLisandro Dalcin /*@
21920f4b53cSBarry Smith PetscPartitionerReset - Resets data structures for the `PetscPartitioner`
220abe9303eSLisandro Dalcin
22120f4b53cSBarry Smith Collective
222abe9303eSLisandro Dalcin
223abe9303eSLisandro Dalcin Input Parameter:
22420f4b53cSBarry Smith . part - the `PetscPartitioner` object to reset
225abe9303eSLisandro Dalcin
226abe9303eSLisandro Dalcin Level: developer
227abe9303eSLisandro Dalcin
228db781477SPatrick Sanan .seealso: `PetscPartitionerSetUp()`, `PetscPartitionerDestroy()`
229abe9303eSLisandro Dalcin @*/
PetscPartitionerReset(PetscPartitioner part)230d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerReset(PetscPartitioner part)
231d71ae5a4SJacob Faibussowitsch {
232abe9303eSLisandro Dalcin PetscFunctionBegin;
233abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
234dbbe0bcdSBarry Smith PetscTryTypeMethod(part, reset);
2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
236abe9303eSLisandro Dalcin }
237abe9303eSLisandro Dalcin
238abe9303eSLisandro Dalcin /*@
23920f4b53cSBarry Smith PetscPartitionerDestroy - Destroys a `PetscPartitioner` object
240abe9303eSLisandro Dalcin
24120f4b53cSBarry Smith Collective
242abe9303eSLisandro Dalcin
243abe9303eSLisandro Dalcin Input Parameter:
24420f4b53cSBarry Smith . part - the `PetscPartitioner` object to destroy
245abe9303eSLisandro Dalcin
246abe9303eSLisandro Dalcin Level: developer
247abe9303eSLisandro Dalcin
248db781477SPatrick Sanan .seealso: `PetscPartitionerView()`
249abe9303eSLisandro Dalcin @*/
PetscPartitionerDestroy(PetscPartitioner * part)250d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
251d71ae5a4SJacob Faibussowitsch {
252abe9303eSLisandro Dalcin PetscFunctionBegin;
2533ba16761SJacob Faibussowitsch if (!*part) PetscFunctionReturn(PETSC_SUCCESS);
254f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*part, PETSCPARTITIONER_CLASSID, 1);
255abe9303eSLisandro Dalcin
256f4f49eeaSPierre Jolivet if (--((PetscObject)*part)->refct > 0) {
2579371c9d4SSatish Balay *part = NULL;
2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
2599371c9d4SSatish Balay }
260f4f49eeaSPierre Jolivet ((PetscObject)*part)->refct = 0;
261abe9303eSLisandro Dalcin
2629566063dSJacob Faibussowitsch PetscCall(PetscPartitionerReset(*part));
263abe9303eSLisandro Dalcin
264648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&(*part)->viewer));
265648c30bcSBarry Smith PetscCall(PetscViewerDestroy(&(*part)->viewerGraph));
266f4f49eeaSPierre Jolivet PetscTryTypeMethod(*part, destroy);
2679566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(part));
2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
269abe9303eSLisandro Dalcin }
270abe9303eSLisandro Dalcin
271abe9303eSLisandro Dalcin /*@
272abe9303eSLisandro Dalcin PetscPartitionerPartition - Partition a graph
273abe9303eSLisandro Dalcin
27420f4b53cSBarry Smith Collective
275abe9303eSLisandro Dalcin
276abe9303eSLisandro Dalcin Input Parameters:
27720f4b53cSBarry Smith + part - The `PetscPartitioner`
278abe9303eSLisandro Dalcin . nparts - Number of partitions
279abe9303eSLisandro Dalcin . numVertices - Number of vertices in the local part of the graph
280abe9303eSLisandro Dalcin . start - row pointers for the local part of the graph (CSR style)
281abe9303eSLisandro Dalcin . adjacency - adjacency list (CSR style)
282dc9a610eSPierre Jolivet . vertexSection - PetscSection describing the absolute weight of each local vertex (can be `NULL`)
28321c92275SMatthew G. Knepley . edgeSection - PetscSection describing the absolute weight of each local edge (can be `NULL`)
284dc9a610eSPierre Jolivet - targetSection - PetscSection describing the absolute weight of each partition (can be `NULL`)
285abe9303eSLisandro Dalcin
286abe9303eSLisandro Dalcin Output Parameters:
28720f4b53cSBarry Smith + partSection - The `PetscSection` giving the division of points by partition
288abe9303eSLisandro Dalcin - partition - The list of points by partition
289abe9303eSLisandro Dalcin
290*3b9d9b65SStefano Zampini Options Database Keys:
291abe9303eSLisandro Dalcin + -petscpartitioner_view - View the partitioner information
292abe9303eSLisandro Dalcin - -petscpartitioner_view_graph - View the graph we are partitioning
293abe9303eSLisandro Dalcin
29420f4b53cSBarry Smith Level: developer
29520f4b53cSBarry Smith
296abe9303eSLisandro Dalcin Notes:
297abe9303eSLisandro Dalcin The chart of the vertexSection (if present) must contain [0,numVertices), with the number of dofs in the section specifying the absolute weight for each vertex.
298abe9303eSLisandro Dalcin The chart of the targetSection (if present) must contain [0,nparts), with the number of dofs in the section specifying the absolute weight for each partition. This information must be the same across processes, PETSc does not check it.
299abe9303eSLisandro Dalcin
300a4e35b19SJacob Faibussowitsch .seealso: `PetscPartitionerCreate()`, `PetscPartitionerSetType()`, `PetscSectionCreate()`, `PetscSectionSetChart()`, `PetscSectionSetDof()`
301abe9303eSLisandro Dalcin @*/
PetscPartitionerPartition(PetscPartitioner part,PetscInt nparts,PetscInt numVertices,PetscInt start[],PetscInt adjacency[],PetscSection vertexSection,PetscSection edgeSection,PetscSection targetSection,PetscSection partSection,IS * partition)30221c92275SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertexSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
303d71ae5a4SJacob Faibussowitsch {
304abe9303eSLisandro Dalcin PetscFunctionBegin;
305abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
306abe9303eSLisandro Dalcin PetscValidLogicalCollectiveInt(part, nparts, 2);
30708401ef6SPierre Jolivet PetscCheck(nparts > 0, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_OUTOFRANGE, "Number of parts must be positive");
30808401ef6SPierre Jolivet PetscCheck(numVertices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices must be non-negative");
309abe9303eSLisandro Dalcin if (numVertices && !part->noGraph) {
3104f572ea9SToby Isaac PetscAssertPointer(start, 4);
3114f572ea9SToby Isaac PetscAssertPointer(start + numVertices, 4);
3124f572ea9SToby Isaac if (start[numVertices]) PetscAssertPointer(adjacency, 5);
313abe9303eSLisandro Dalcin }
314abe9303eSLisandro Dalcin if (vertexSection) {
315abe9303eSLisandro Dalcin PetscInt s, e;
316abe9303eSLisandro Dalcin
317abe9303eSLisandro Dalcin PetscValidHeaderSpecific(vertexSection, PETSC_SECTION_CLASSID, 6);
3189566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(vertexSection, &s, &e));
3191dca8a05SBarry Smith PetscCheck(s <= 0 && e >= numVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid vertexSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
320abe9303eSLisandro Dalcin }
32121c92275SMatthew G. Knepley if (edgeSection) {
32221c92275SMatthew G. Knepley PetscInt s, e;
32321c92275SMatthew G. Knepley
32421c92275SMatthew G. Knepley PetscValidHeaderSpecific(edgeSection, PETSC_SECTION_CLASSID, 7);
32521c92275SMatthew G. Knepley PetscCall(PetscSectionGetChart(edgeSection, &s, &e));
32621c92275SMatthew G. Knepley PetscCheck(s <= 0 && e >= start[numVertices], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid edgeSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
32721c92275SMatthew G. Knepley }
328abe9303eSLisandro Dalcin if (targetSection) {
329abe9303eSLisandro Dalcin PetscInt s, e;
330abe9303eSLisandro Dalcin
33121c92275SMatthew G. Knepley PetscValidHeaderSpecific(targetSection, PETSC_SECTION_CLASSID, 8);
3329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(targetSection, &s, &e));
3331dca8a05SBarry Smith PetscCheck(s <= 0 && e >= nparts, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid targetSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
334abe9303eSLisandro Dalcin }
33521c92275SMatthew G. Knepley PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 9);
33621c92275SMatthew G. Knepley PetscAssertPointer(partition, 10);
337abe9303eSLisandro Dalcin
3389566063dSJacob Faibussowitsch PetscCall(PetscSectionReset(partSection));
3399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(partSection, 0, nparts));
340abe9303eSLisandro Dalcin if (nparts == 1) { /* quick */
3419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(partSection, 0, numVertices));
3429566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)part), numVertices, 0, 1, partition));
34321c92275SMatthew G. Knepley } else PetscUseTypeMethod(part, partition, nparts, numVertices, start, adjacency, vertexSection, edgeSection, targetSection, partSection, partition);
3449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(partSection));
345abe9303eSLisandro Dalcin if (part->viewerGraph) {
346abe9303eSLisandro Dalcin PetscViewer viewer = part->viewerGraph;
347abe9303eSLisandro Dalcin PetscBool isascii;
348abe9303eSLisandro Dalcin PetscInt v, i;
349abe9303eSLisandro Dalcin PetscMPIInt rank;
350abe9303eSLisandro Dalcin
3519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
3529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
353abe9303eSLisandro Dalcin if (isascii) {
3549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer));
35563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]Nv: %" PetscInt_FMT "\n", rank, numVertices));
356abe9303eSLisandro Dalcin for (v = 0; v < numVertices; ++v) {
357abe9303eSLisandro Dalcin const PetscInt s = start[v];
358abe9303eSLisandro Dalcin const PetscInt e = start[v + 1];
359abe9303eSLisandro Dalcin
3609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] ", rank));
36163a3b9bcSJacob Faibussowitsch for (i = s; i < e; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT " ", adjacency[i]));
36263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%" PetscInt_FMT "-%" PetscInt_FMT ")\n", s, e));
363abe9303eSLisandro Dalcin }
3649566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer));
3659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer));
366abe9303eSLisandro Dalcin }
367abe9303eSLisandro Dalcin }
368*3b9d9b65SStefano Zampini if (part->viewer) {
369*3b9d9b65SStefano Zampini PetscCall(PetscViewerPushFormat(part->viewer, part->viewerFmt));
370*3b9d9b65SStefano Zampini PetscCall(PetscPartitionerView(part, part->viewer));
371*3b9d9b65SStefano Zampini PetscCall(PetscViewerPopFormat(part->viewer));
372*3b9d9b65SStefano Zampini }
3733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
374abe9303eSLisandro Dalcin }
375abe9303eSLisandro Dalcin
376abe9303eSLisandro Dalcin /*@
37720f4b53cSBarry Smith PetscPartitionerCreate - Creates an empty `PetscPartitioner` object. The type can then be set with `PetscPartitionerSetType()`.
378abe9303eSLisandro Dalcin
379abe9303eSLisandro Dalcin Collective
380abe9303eSLisandro Dalcin
381abe9303eSLisandro Dalcin Input Parameter:
38220f4b53cSBarry Smith . comm - The communicator for the `PetscPartitioner` object
383abe9303eSLisandro Dalcin
384abe9303eSLisandro Dalcin Output Parameter:
38520f4b53cSBarry Smith . part - The `PetscPartitioner` object
386abe9303eSLisandro Dalcin
387abe9303eSLisandro Dalcin Level: beginner
388abe9303eSLisandro Dalcin
389db781477SPatrick Sanan .seealso: `PetscPartitionerSetType()`, `PetscPartitionerDestroy()`
390abe9303eSLisandro Dalcin @*/
PetscPartitionerCreate(MPI_Comm comm,PetscPartitioner * part)391d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
392d71ae5a4SJacob Faibussowitsch {
393abe9303eSLisandro Dalcin PetscPartitioner p;
394abe9303eSLisandro Dalcin const char *partitionerType = NULL;
395abe9303eSLisandro Dalcin
396abe9303eSLisandro Dalcin PetscFunctionBegin;
3974f572ea9SToby Isaac PetscAssertPointer(part, 2);
398abe9303eSLisandro Dalcin *part = NULL;
3999566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitializePackage());
400abe9303eSLisandro Dalcin
4019566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView));
4029566063dSJacob Faibussowitsch PetscCall(PetscPartitionerGetDefaultType(comm, &partitionerType));
4039566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(p, partitionerType));
404abe9303eSLisandro Dalcin
405abe9303eSLisandro Dalcin p->usevwgt = PETSC_TRUE;
406*3b9d9b65SStefano Zampini p->printHeader = PETSC_TRUE;
407abe9303eSLisandro Dalcin
408abe9303eSLisandro Dalcin *part = p;
4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
410abe9303eSLisandro Dalcin }
411