xref: /petsc/src/dm/partitioner/impls/chaco/partchaco.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
2abe9303eSLisandro Dalcin 
3abe9303eSLisandro Dalcin PetscBool  ChacoPartitionerCite       = PETSC_FALSE;
49371c9d4SSatish Balay const char ChacoPartitionerCitation[] = "@inproceedings{Chaco95,\n"
5abe9303eSLisandro Dalcin                                         "  author    = {Bruce Hendrickson and Robert Leland},\n"
6abe9303eSLisandro Dalcin                                         "  title     = {A multilevel algorithm for partitioning graphs},\n"
7abe9303eSLisandro Dalcin                                         "  booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)},"
8abe9303eSLisandro Dalcin                                         "  isbn      = {0-89791-816-9},\n"
9abe9303eSLisandro Dalcin                                         "  pages     = {28},\n"
10abe9303eSLisandro Dalcin                                         "  doi       = {https://doi.acm.org/10.1145/224170.224228},\n"
11abe9303eSLisandro Dalcin                                         "  publisher = {ACM Press},\n"
12abe9303eSLisandro Dalcin                                         "  address   = {New York},\n"
13abe9303eSLisandro Dalcin                                         "  year      = {1995}\n"
14abe9303eSLisandro Dalcin                                         "}\n";
15abe9303eSLisandro Dalcin 
16abe9303eSLisandro Dalcin typedef struct {
17abe9303eSLisandro Dalcin   PetscInt dummy;
18abe9303eSLisandro Dalcin } PetscPartitioner_Chaco;
19abe9303eSLisandro Dalcin 
209371c9d4SSatish Balay static PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part) {
21abe9303eSLisandro Dalcin   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *)part->data;
22abe9303eSLisandro Dalcin 
23abe9303eSLisandro Dalcin   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(PetscFree(p));
25abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
26abe9303eSLisandro Dalcin }
27abe9303eSLisandro Dalcin 
289371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Chaco_ASCII(PetscPartitioner part, PetscViewer viewer) {
29abe9303eSLisandro Dalcin   PetscFunctionBegin;
30abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
31abe9303eSLisandro Dalcin }
32abe9303eSLisandro Dalcin 
339371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer) {
34abe9303eSLisandro Dalcin   PetscBool iascii;
35abe9303eSLisandro Dalcin 
36abe9303eSLisandro Dalcin   PetscFunctionBegin;
37abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
38abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
399566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
409566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Chaco_ASCII(part, viewer));
41abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
42abe9303eSLisandro Dalcin }
43abe9303eSLisandro Dalcin 
44abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
45abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_UNISTD_H)
46abe9303eSLisandro Dalcin #include <unistd.h>
47abe9303eSLisandro Dalcin #endif
48abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
49abe9303eSLisandro Dalcin #include <chaco.h>
50abe9303eSLisandro Dalcin #else
51abe9303eSLisandro Dalcin /* Older versions of Chaco do not have an include file */
529371c9d4SSatish Balay PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts, float *ewgts, float *x, float *y, float *z, char *outassignname, char *outfilename, short *assignment, int architecture, int ndims_tot, int mesh_dims[3], double *goal, int global_method, int local_method, int rqi_flag, int vmax, int ndims, double eigtol, long seed);
53abe9303eSLisandro Dalcin #endif
54abe9303eSLisandro Dalcin extern int FREE_GRAPH;
55abe9303eSLisandro Dalcin #endif
56abe9303eSLisandro Dalcin 
579371c9d4SSatish Balay static PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) {
58abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
599371c9d4SSatish Balay   enum {
609371c9d4SSatish Balay     DEFAULT_METHOD  = 1,
619371c9d4SSatish Balay     INERTIAL_METHOD = 3
629371c9d4SSatish Balay   };
63abe9303eSLisandro Dalcin   MPI_Comm comm;
64abe9303eSLisandro Dalcin   int      nvtxs = numVertices;            /* number of vertices in full graph */
65abe9303eSLisandro Dalcin   int     *vwgts = NULL;                   /* weights for all vertices */
66abe9303eSLisandro Dalcin   float   *ewgts = NULL;                   /* weights for all edges */
67abe9303eSLisandro Dalcin   float   *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
68abe9303eSLisandro Dalcin   char    *outassignname = NULL;           /*  name of assignment output file */
69abe9303eSLisandro Dalcin   char    *outfilename   = NULL;           /* output file name */
70abe9303eSLisandro Dalcin   int      architecture  = 1;              /* 0 => hypercube, d => d-dimensional mesh */
71abe9303eSLisandro Dalcin   int      ndims_tot     = 0;              /* total number of cube dimensions to divide */
72abe9303eSLisandro Dalcin   int      mesh_dims[3];                   /* dimensions of mesh of processors */
73abe9303eSLisandro Dalcin   double  *goal          = NULL;           /* desired set sizes for each set */
74abe9303eSLisandro Dalcin   int      global_method = 1;              /* global partitioning algorithm */
75abe9303eSLisandro Dalcin   int      local_method  = 1;              /* local partitioning algorithm */
76abe9303eSLisandro Dalcin   int      rqi_flag      = 0;              /* should I use RQI/Symmlq eigensolver? */
77abe9303eSLisandro Dalcin   int      vmax          = 200;            /* how many vertices to coarsen down to? */
78abe9303eSLisandro Dalcin   int      ndims         = 1;              /* number of eigenvectors (2^d sets) */
79abe9303eSLisandro Dalcin   double   eigtol        = 0.001;          /* tolerance on eigenvectors */
80abe9303eSLisandro Dalcin   long     seed          = 123636512;      /* for random graph mutations */
81abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
82abe9303eSLisandro Dalcin   int *assignment; /* Output partition */
83abe9303eSLisandro Dalcin #else
84abe9303eSLisandro Dalcin   short int *assignment; /* Output partition */
85abe9303eSLisandro Dalcin #endif
86abe9303eSLisandro Dalcin   int       fd_stdout, fd_pipe[2];
87abe9303eSLisandro Dalcin   PetscInt *points;
88abe9303eSLisandro Dalcin   int       i, v, p;
89d0609cedSBarry Smith   int       err;
90abe9303eSLisandro Dalcin 
91abe9303eSLisandro Dalcin   PetscFunctionBegin;
929566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
93abe9303eSLisandro Dalcin   if (PetscDefined(USE_DEBUG)) {
94abe9303eSLisandro Dalcin     int       ival, isum;
95abe9303eSLisandro Dalcin     PetscBool distributed;
96abe9303eSLisandro Dalcin 
97abe9303eSLisandro Dalcin     ival = (numVertices > 0);
989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Allreduce(&ival, &isum, 1, MPI_INT, MPI_SUM, comm));
99abe9303eSLisandro Dalcin     distributed = (isum > 1) ? PETSC_TRUE : PETSC_FALSE;
10028b400f6SJacob Faibussowitsch     PetscCheck(!distributed, comm, PETSC_ERR_SUP, "Chaco cannot partition a distributed graph");
101abe9303eSLisandro Dalcin   }
102abe9303eSLisandro Dalcin   if (!numVertices) { /* distributed case, return if not holding the graph */
1039566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition));
104abe9303eSLisandro Dalcin     PetscFunctionReturn(0);
105abe9303eSLisandro Dalcin   }
106abe9303eSLisandro Dalcin   FREE_GRAPH = 0; /* Do not let Chaco free my memory */
107abe9303eSLisandro Dalcin   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
108abe9303eSLisandro Dalcin 
109e08b1d6dSBarry Smith   /* code would use manager.createCellCoordinates(nvtxs, &x, &y, &z); */
110e08b1d6dSBarry Smith   PetscCheck(global_method != INERTIAL_METHOD, PETSC_COMM_SELF, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
111e08b1d6dSBarry Smith 
112abe9303eSLisandro Dalcin   mesh_dims[0] = nparts;
113abe9303eSLisandro Dalcin   mesh_dims[1] = 1;
114abe9303eSLisandro Dalcin   mesh_dims[2] = 1;
1159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nvtxs, &assignment));
116abe9303eSLisandro Dalcin   /* Chaco outputs to stdout. We redirect this to a buffer. */
117abe9303eSLisandro Dalcin   /* TODO: check error codes for UNIX calls */
118abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_UNISTD_H)
119abe9303eSLisandro Dalcin   {
120abe9303eSLisandro Dalcin     int piperet;
121abe9303eSLisandro Dalcin     piperet = pipe(fd_pipe);
12228b400f6SJacob Faibussowitsch     PetscCheck(!piperet, PETSC_COMM_SELF, PETSC_ERR_SYS, "Could not create pipe");
123abe9303eSLisandro Dalcin     fd_stdout = dup(1);
124abe9303eSLisandro Dalcin     close(1);
125abe9303eSLisandro Dalcin     dup2(fd_pipe[1], 1);
126abe9303eSLisandro Dalcin   }
127abe9303eSLisandro Dalcin #endif
1289566063dSJacob Faibussowitsch   if (part->usevwgt) PetscCall(PetscInfo(part, "PETSCPARTITIONERCHACO ignores vertex weights\n"));
1299371c9d4SSatish Balay   err = interface(nvtxs, (int *)start, (int *)adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename, assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims, eigtol, seed);
130abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_UNISTD_H)
131abe9303eSLisandro Dalcin   {
132abe9303eSLisandro Dalcin     char msgLog[10000];
133abe9303eSLisandro Dalcin     int  count;
134abe9303eSLisandro Dalcin 
135abe9303eSLisandro Dalcin     fflush(stdout);
136abe9303eSLisandro Dalcin     count = read(fd_pipe[0], msgLog, (10000 - 1) * sizeof(char));
137abe9303eSLisandro Dalcin     if (count < 0) count = 0;
138abe9303eSLisandro Dalcin     msgLog[count] = 0;
139abe9303eSLisandro Dalcin     close(1);
140abe9303eSLisandro Dalcin     dup2(fd_stdout, 1);
141abe9303eSLisandro Dalcin     close(fd_stdout);
142abe9303eSLisandro Dalcin     close(fd_pipe[0]);
143abe9303eSLisandro Dalcin     close(fd_pipe[1]);
144d0609cedSBarry Smith     PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
145abe9303eSLisandro Dalcin   }
146abe9303eSLisandro Dalcin #else
147d0609cedSBarry Smith   PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", "error in stdout");
148abe9303eSLisandro Dalcin #endif
149abe9303eSLisandro Dalcin   /* Convert to PetscSection+IS */
15048a46eb9SPierre Jolivet   for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1));
1519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nvtxs, &points));
152abe9303eSLisandro Dalcin   for (p = 0, i = 0; p < nparts; ++p) {
153abe9303eSLisandro Dalcin     for (v = 0; v < nvtxs; ++v) {
154abe9303eSLisandro Dalcin       if (assignment[v] == p) points[i++] = v;
155abe9303eSLisandro Dalcin     }
156abe9303eSLisandro Dalcin   }
15763a3b9bcSJacob Faibussowitsch   PetscCheck(i == nvtxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs);
1589566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition));
159e08b1d6dSBarry Smith   /* code would use manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
160e08b1d6dSBarry Smith   PetscCheck(global_method != INERTIAL_METHOD, PETSC_COMM_SELF, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
1619566063dSJacob Faibussowitsch   PetscCall(PetscFree(assignment));
162abe9303eSLisandro Dalcin   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
163abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
164abe9303eSLisandro Dalcin #else
165abe9303eSLisandro Dalcin   SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
166abe9303eSLisandro Dalcin #endif
167abe9303eSLisandro Dalcin }
168abe9303eSLisandro Dalcin 
1699371c9d4SSatish Balay static PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part) {
170abe9303eSLisandro Dalcin   PetscFunctionBegin;
171abe9303eSLisandro Dalcin   part->noGraph        = PETSC_FALSE;
172abe9303eSLisandro Dalcin   part->ops->view      = PetscPartitionerView_Chaco;
173abe9303eSLisandro Dalcin   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
174abe9303eSLisandro Dalcin   part->ops->partition = PetscPartitionerPartition_Chaco;
175abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
176abe9303eSLisandro Dalcin }
177abe9303eSLisandro Dalcin 
178abe9303eSLisandro Dalcin /*MC
179abe9303eSLisandro Dalcin   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
180abe9303eSLisandro Dalcin 
181abe9303eSLisandro Dalcin   Level: intermediate
182abe9303eSLisandro Dalcin 
183db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
184abe9303eSLisandro Dalcin M*/
185abe9303eSLisandro Dalcin 
1869371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part) {
187abe9303eSLisandro Dalcin   PetscPartitioner_Chaco *p;
188abe9303eSLisandro Dalcin 
189abe9303eSLisandro Dalcin   PetscFunctionBegin;
190abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
191*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
192abe9303eSLisandro Dalcin   part->data = p;
193abe9303eSLisandro Dalcin 
1949566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Chaco(part));
1959566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionerCite));
196abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
197abe9303eSLisandro Dalcin }
198