Actual source code: partitionerimpl.h
1: #pragma once
3: #include <petscviewer.h>
4: #include <petscpartitioner.h>
5: #include <petsc/private/petscimpl.h>
7: PETSC_EXTERN PetscBool PetscPartitionerRegisterAllCalled;
8: PETSC_EXTERN PetscErrorCode PetscPartitionerRegisterAll(void);
10: typedef struct _PetscPartitionerOps *PetscPartitionerOps;
11: struct _PetscPartitionerOps {
12: PetscErrorCode (*setfromoptions)(PetscPartitioner, PetscOptionItems);
13: PetscErrorCode (*setup)(PetscPartitioner);
14: PetscErrorCode (*reset)(PetscPartitioner);
15: PetscErrorCode (*view)(PetscPartitioner, PetscViewer);
16: PetscErrorCode (*destroy)(PetscPartitioner);
17: PetscErrorCode (*partition)(PetscPartitioner, PetscInt, PetscInt, PetscInt[], PetscInt[], PetscSection, PetscSection, PetscSection, PetscSection, IS *);
18: };
20: struct _p_PetscPartitioner {
21: PETSCHEADER(struct _PetscPartitionerOps);
22: void *data; /* Implementation object */
23: PetscInt height; /* Height of points to partition into non-overlapping subsets */
24: PetscInt edgeCut; /* The number of edge cut by the partition */
25: PetscReal balance; /* The maximum partition size divided by the minimum size */
27: PetscBool printHeader;
28: PetscViewer viewer, viewerGraph;
29: PetscViewerFormat viewerFmt;
31: PetscBool noGraph; /* if true, the partitioner does not need the connectivity graph, only the number of local vertices */
32: PetscBool usevwgt; /* if true, the partitioner looks at the local section vertSection to weight the vertices of the graph */
33: PetscBool useewgt; /* if true, the partitioner looks at the topology to weight the edges of the graph */
34: };
36: /* All levels > 8 logged in 8-th level */
37: #define PETSCPARTITIONER_MS_MAXSTAGE 8
38: #define PETSCPARTITIONER_MS_NUMSTAGE (PETSCPARTITIONER_MS_MAXSTAGE + 1)
39: PETSC_EXTERN PetscLogEvent PetscPartitioner_MS_SetUp;
40: PETSC_EXTERN PetscLogEvent PetscPartitioner_MS_Stage[PETSCPARTITIONER_MS_NUMSTAGE];
41: PETSC_EXTERN PetscErrorCode PetscPartitionerMultistageGetStages_Multistage(PetscPartitioner, PetscInt *, MPI_Group *[]);
42: PETSC_EXTERN PetscErrorCode PetscPartitionerMultistageSetStage_Multistage(PetscPartitioner, PetscInt, PetscObject);
43: PETSC_EXTERN PetscErrorCode PetscPartitionerMultistageGetStage_Multistage(PetscPartitioner, PetscInt *, PetscObject *);