Actual source code: partitioner.c

  1: #include <petsc/private/partitionerimpl.h>

  3: /*@C
  4:   PetscPartitionerSetType - Builds a particular `PetscPartitioner`

  6:   Collective

  8:   Input Parameters:
  9: + part - The `PetscPartitioner` object
 10: - name - The kind of partitioner

 12:   Options Database Key:
 13: . -petscpartitioner_type <type> - Sets the `PetscPartitioner` type

 15:   Level: intermediate

 17:   Note:
 18: .vb
 19:  PETSCPARTITIONERCHACO    - The Chaco partitioner (--download-chaco)
 20:  PETSCPARTITIONERPARMETIS - The ParMetis partitioner (--download-parmetis)
 21:  PETSCPARTITIONERSHELL    - A shell partitioner implemented by the user
 22:  PETSCPARTITIONERSIMPLE   - A simple partitioner that divides cells into equal, contiguous chunks
 23:  PETSCPARTITIONERGATHER   - Gathers all cells onto process 0
 24: .ve

 26: .seealso: `PetscPartitionerGetType()`, `PetscPartitionerCreate()`
 27: @*/
 28: PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
 29: {
 30:   PetscErrorCode (*r)(PetscPartitioner);
 31:   PetscBool match;

 33:   PetscFunctionBegin;
 35:   PetscCall(PetscObjectTypeCompare((PetscObject)part, name, &match));
 36:   if (match) PetscFunctionReturn(PETSC_SUCCESS);

 38:   PetscCall(PetscPartitionerRegisterAll());
 39:   PetscCall(PetscFunctionListFind(PetscPartitionerList, name, &r));
 40:   PetscCheck(r, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);

 42:   PetscTryTypeMethod(part, destroy);
 43:   part->noGraph = PETSC_FALSE;
 44:   PetscCall(PetscMemzero(part->ops, sizeof(*part->ops)));
 45:   PetscCall(PetscObjectChangeTypeName((PetscObject)part, name));
 46:   PetscCall((*r)(part));
 47:   PetscFunctionReturn(PETSC_SUCCESS);
 48: }

 50: /*@C
 51:   PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.

 53:   Not Collective

 55:   Input Parameter:
 56: . part - The PetscPartitioner

 58:   Output Parameter:
 59: . name - The PetscPartitioner type name

 61:   Level: intermediate

 63: .seealso: `PetscPartitionerSetType()`, `PetscPartitionerCreate()`
 64: @*/
 65: PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
 66: {
 67:   PetscFunctionBegin;
 69:   PetscAssertPointer(name, 2);
 70:   *name = ((PetscObject)part)->type_name;
 71:   PetscFunctionReturn(PETSC_SUCCESS);
 72: }

 74: /*@C
 75:   PetscPartitionerViewFromOptions - View a `PetscPartitioner` object based on options in the options database

 77:   Collective

 79:   Input Parameters:
 80: + A    - the `PetscPartitioner` object
 81: . obj  - Optional `PetscObject` that provides the options prefix
 82: - name - command line option

 84:   Level: intermediate

 86:   Note:
 87:   See `PetscObjectViewFromOptions()` for the various forms of viewers that may be used

 89: .seealso: `PetscPartitionerView()`, `PetscObjectViewFromOptions()`
 90: @*/
 91: PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner A, PetscObject obj, const char name[])
 92: {
 93:   PetscFunctionBegin;
 95:   PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
 96:   PetscFunctionReturn(PETSC_SUCCESS);
 97: }

 99: /*@
100:   PetscPartitionerView - Views a `PetscPartitioner`

102:   Collective

104:   Input Parameters:
105: + part - the `PetscPartitioner` object to view
106: - v    - the viewer

108:   Level: developer

110: .seealso: `PetscPartitionerDestroy()`
111: @*/
112: PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
113: {
114:   PetscMPIInt size;
115:   PetscBool   isascii;

117:   PetscFunctionBegin;
119:   if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)part), &v));
120:   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &isascii));
121:   if (isascii) {
122:     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)part), &size));
123:     PetscCall(PetscViewerASCIIPrintf(v, "Graph Partitioner: %d MPI Process%s\n", size, size > 1 ? "es" : ""));
124:     PetscCall(PetscViewerASCIIPrintf(v, "  type: %s\n", ((PetscObject)part)->type_name));
125:     PetscCall(PetscViewerASCIIPrintf(v, "  edge cut: %" PetscInt_FMT "\n", part->edgeCut));
126:     PetscCall(PetscViewerASCIIPrintf(v, "  balance: %.2g\n", (double)part->balance));
127:     PetscCall(PetscViewerASCIIPrintf(v, "  use vertex weights: %d\n", part->usevwgt));
128:     PetscCall(PetscViewerASCIIPrintf(v, "  use edge weights: %d\n", part->useewgt));
129:   }
130:   PetscTryTypeMethod(part, view, v);
131:   PetscFunctionReturn(PETSC_SUCCESS);
132: }

134: static PetscErrorCode PetscPartitionerGetDefaultType(MPI_Comm comm, const char **defaultType)
135: {
136:   PetscMPIInt size;

138:   PetscFunctionBegin;
139:   PetscCallMPI(MPI_Comm_size(comm, &size));
140:   if (size == 1) {
141:     *defaultType = PETSCPARTITIONERSIMPLE;
142:   } else {
143: #if defined(PETSC_HAVE_PARMETIS)
144:     *defaultType = PETSCPARTITIONERPARMETIS;
145: #elif defined(PETSC_HAVE_PTSCOTCH)
146:     *defaultType = PETSCPARTITIONERPTSCOTCH;
147: #elif defined(PETSC_HAVE_CHACO)
148:     *defaultType = PETSCPARTITIONERCHACO;
149: #else
150:     *defaultType = PETSCPARTITIONERSIMPLE;
151: #endif
152:   }
153:   PetscFunctionReturn(PETSC_SUCCESS);
154: }

156: /*@
157:   PetscPartitionerSetFromOptions - sets parameters in a `PetscPartitioner` from the options database

159:   Collective

161:   Input Parameter:
162: . part - the `PetscPartitioner` object to set options for

164:   Options Database Keys:
165: + -petscpartitioner_type <type>        - Sets the `PetscPartitioner` type; use -help for a list of available types
166: . -petscpartitioner_use_vertex_weights - Uses weights associated with the graph vertices
167: - -petscpartitioner_view_graph         - View the graph each time PetscPartitionerPartition is called. Viewer can be customized, see `PetscOptionsGetViewer()`

169:   Level: developer

171: .seealso: `PetscPartitionerView()`, `PetscPartitionerSetType()`, `PetscPartitionerPartition()`
172: @*/
173: PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
174: {
175:   const char *currentType = NULL;
176:   char        name[256];
177:   PetscBool   flg;

179:   PetscFunctionBegin;
181:   PetscObjectOptionsBegin((PetscObject)part);
182:   PetscCall(PetscPartitionerGetType(part, &currentType));
183:   PetscCall(PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, currentType, name, sizeof(name), &flg));
184:   if (flg) PetscCall(PetscPartitionerSetType(part, name));
185:   PetscCall(PetscOptionsBool("-petscpartitioner_use_vertex_weights", "Use vertex weights", "", part->usevwgt, &part->usevwgt, NULL));
186:   PetscCall(PetscOptionsBool("-petscpartitioner_use_edge_weights", "Use edge weights", "", part->useewgt, &part->useewgt, NULL));
187:   PetscTryTypeMethod(part, setfromoptions, PetscOptionsObject);
188:   PetscCall(PetscOptionsRestoreViewer(&part->viewer));
189:   PetscCall(PetscOptionsRestoreViewer(&part->viewerGraph));
190:   PetscCall(PetscOptionsGetViewer(((PetscObject)part)->comm, ((PetscObject)part)->options, ((PetscObject)part)->prefix, "-petscpartitioner_view", &part->viewer, NULL, NULL));
191:   PetscCall(PetscOptionsGetViewer(((PetscObject)part)->comm, ((PetscObject)part)->options, ((PetscObject)part)->prefix, "-petscpartitioner_view_graph", &part->viewerGraph, NULL, &part->viewGraph));
192:   /* process any options handlers added with PetscObjectAddOptionsHandler() */
193:   PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)part, PetscOptionsObject));
194:   PetscOptionsEnd();
195:   PetscFunctionReturn(PETSC_SUCCESS);
196: }

198: /*@
199:   PetscPartitionerSetUp - Construct data structures for the `PetscPartitioner`

201:   Collective

203:   Input Parameter:
204: . part - the `PetscPartitioner` object to setup

206:   Level: developer

208: .seealso: `PetscPartitionerView()`, `PetscPartitionerDestroy()`
209: @*/
210: PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
211: {
212:   PetscFunctionBegin;
214:   PetscTryTypeMethod(part, setup);
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: /*@
219:   PetscPartitionerReset - Resets data structures for the `PetscPartitioner`

221:   Collective

223:   Input Parameter:
224: . part - the `PetscPartitioner` object to reset

226:   Level: developer

228: .seealso: `PetscPartitionerSetUp()`, `PetscPartitionerDestroy()`
229: @*/
230: PetscErrorCode PetscPartitionerReset(PetscPartitioner part)
231: {
232:   PetscFunctionBegin;
234:   PetscTryTypeMethod(part, reset);
235:   PetscFunctionReturn(PETSC_SUCCESS);
236: }

238: /*@
239:   PetscPartitionerDestroy - Destroys a `PetscPartitioner` object

241:   Collective

243:   Input Parameter:
244: . part - the `PetscPartitioner` object to destroy

246:   Level: developer

248: .seealso: `PetscPartitionerView()`
249: @*/
250: PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
251: {
252:   PetscFunctionBegin;
253:   if (!*part) PetscFunctionReturn(PETSC_SUCCESS);

256:   if (--((PetscObject)*part)->refct > 0) {
257:     *part = NULL;
258:     PetscFunctionReturn(PETSC_SUCCESS);
259:   }
260:   ((PetscObject)*part)->refct = 0;

262:   PetscCall(PetscPartitionerReset(*part));

264:   PetscCall(PetscOptionsRestoreViewer(&(*part)->viewer));
265:   PetscCall(PetscOptionsRestoreViewer(&(*part)->viewerGraph));
266:   PetscTryTypeMethod(*part, destroy);
267:   PetscCall(PetscHeaderDestroy(part));
268:   PetscFunctionReturn(PETSC_SUCCESS);
269: }

271: /*@
272:   PetscPartitionerPartition - Partition a graph

274:   Collective

276:   Input Parameters:
277: + part          - The `PetscPartitioner`
278: . nparts        - Number of partitions
279: . numVertices   - Number of vertices in the local part of the graph
280: . start         - row pointers for the local part of the graph (CSR style)
281: . adjacency     - adjacency list (CSR style)
282: . vertexSection - PetscSection describing the absolute weight of each local vertex (can be `NULL`)
283: . edgeSection   - PetscSection describing the absolute weight of each local edge (can be `NULL`)
284: - targetSection - PetscSection describing the absolute weight of each partition (can be `NULL`)

286:   Output Parameters:
287: + partSection - The `PetscSection` giving the division of points by partition
288: - partition   - The list of points by partition

290:   Options Databasen Keys:
291: + -petscpartitioner_view       - View the partitioner information
292: - -petscpartitioner_view_graph - View the graph we are partitioning

294:   Level: developer

296:   Notes:
297:   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.
298:   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.

300: .seealso: `PetscPartitionerCreate()`, `PetscPartitionerSetType()`, `PetscSectionCreate()`, `PetscSectionSetChart()`, `PetscSectionSetDof()`
301: @*/
302: PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertexSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
303: {
304:   PetscFunctionBegin;
307:   PetscCheck(nparts > 0, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_OUTOFRANGE, "Number of parts must be positive");
308:   PetscCheck(numVertices >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices must be non-negative");
309:   if (numVertices && !part->noGraph) {
310:     PetscAssertPointer(start, 4);
311:     PetscAssertPointer(start + numVertices, 4);
312:     if (start[numVertices]) PetscAssertPointer(adjacency, 5);
313:   }
314:   if (vertexSection) {
315:     PetscInt s, e;

318:     PetscCall(PetscSectionGetChart(vertexSection, &s, &e));
319:     PetscCheck(s <= 0 && e >= numVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid vertexSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
320:   }
321:   if (edgeSection) {
322:     PetscInt s, e;

325:     PetscCall(PetscSectionGetChart(edgeSection, &s, &e));
326:     PetscCheck(s <= 0 && e >= start[numVertices], PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid edgeSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
327:   }
328:   if (targetSection) {
329:     PetscInt s, e;

332:     PetscCall(PetscSectionGetChart(targetSection, &s, &e));
333:     PetscCheck(s <= 0 && e >= nparts, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid targetSection chart [%" PetscInt_FMT ",%" PetscInt_FMT ")", s, e);
334:   }
336:   PetscAssertPointer(partition, 10);

338:   PetscCall(PetscSectionReset(partSection));
339:   PetscCall(PetscSectionSetChart(partSection, 0, nparts));
340:   if (nparts == 1) { /* quick */
341:     PetscCall(PetscSectionSetDof(partSection, 0, numVertices));
342:     PetscCall(ISCreateStride(PetscObjectComm((PetscObject)part), numVertices, 0, 1, partition));
343:   } else PetscUseTypeMethod(part, partition, nparts, numVertices, start, adjacency, vertexSection, edgeSection, targetSection, partSection, partition);
344:   PetscCall(PetscSectionSetUp(partSection));
345:   if (part->viewerGraph) {
346:     PetscViewer viewer = part->viewerGraph;
347:     PetscBool   isascii;
348:     PetscInt    v, i;
349:     PetscMPIInt rank;

351:     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
352:     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
353:     if (isascii) {
354:       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
355:       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]Nv: %" PetscInt_FMT "\n", rank, numVertices));
356:       for (v = 0; v < numVertices; ++v) {
357:         const PetscInt s = start[v];
358:         const PetscInt e = start[v + 1];

360:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d]  ", rank));
361:         for (i = s; i < e; ++i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT " ", adjacency[i]));
362:         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%" PetscInt_FMT "-%" PetscInt_FMT ")\n", s, e));
363:       }
364:       PetscCall(PetscViewerFlush(viewer));
365:       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
366:     }
367:   }
368:   if (part->viewer) PetscCall(PetscPartitionerView(part, part->viewer));
369:   PetscFunctionReturn(PETSC_SUCCESS);
370: }

372: /*@
373:   PetscPartitionerCreate - Creates an empty `PetscPartitioner` object. The type can then be set with `PetscPartitionerSetType()`.

375:   Collective

377:   Input Parameter:
378: . comm - The communicator for the `PetscPartitioner` object

380:   Output Parameter:
381: . part - The `PetscPartitioner` object

383:   Level: beginner

385: .seealso: `PetscPartitionerSetType()`, `PetscPartitionerDestroy()`
386: @*/
387: PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
388: {
389:   PetscPartitioner p;
390:   const char      *partitionerType = NULL;

392:   PetscFunctionBegin;
393:   PetscAssertPointer(part, 2);
394:   *part = NULL;
395:   PetscCall(PetscPartitionerInitializePackage());

397:   PetscCall(PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView));
398:   PetscCall(PetscPartitionerGetDefaultType(comm, &partitionerType));
399:   PetscCall(PetscPartitionerSetType(p, partitionerType));

401:   p->edgeCut = 0;
402:   p->balance = 0.0;
403:   p->usevwgt = PETSC_TRUE;

405:   *part = p;
406:   PetscFunctionReturn(PETSC_SUCCESS);
407: }