xref: /petsc/src/dm/partitioner/impls/chaco/partchaco.c (revision 6a210b70a61472cf8733db59e8a3fa68f6578c01)
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 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
21d71ae5a4SJacob Faibussowitsch {
22abe9303eSLisandro Dalcin   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *)part->data;
23abe9303eSLisandro Dalcin 
24abe9303eSLisandro Dalcin   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(p));
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27abe9303eSLisandro Dalcin }
28abe9303eSLisandro Dalcin 
29d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Chaco_ASCII(PetscPartitioner part, PetscViewer viewer)
30d71ae5a4SJacob Faibussowitsch {
31abe9303eSLisandro Dalcin   PetscFunctionBegin;
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33abe9303eSLisandro Dalcin }
34abe9303eSLisandro Dalcin 
35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
36d71ae5a4SJacob Faibussowitsch {
37abe9303eSLisandro Dalcin   PetscBool iascii;
38abe9303eSLisandro Dalcin 
39abe9303eSLisandro Dalcin   PetscFunctionBegin;
40abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
41abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
439566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Chaco_ASCII(part, viewer));
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
45abe9303eSLisandro Dalcin }
46abe9303eSLisandro Dalcin 
47abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
48abe9303eSLisandro Dalcin   #if defined(PETSC_HAVE_UNISTD_H)
49abe9303eSLisandro Dalcin     #include <unistd.h>
50abe9303eSLisandro Dalcin   #endif
51abe9303eSLisandro Dalcin   #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
52abe9303eSLisandro Dalcin     #include <chaco.h>
53abe9303eSLisandro Dalcin   #else
54abe9303eSLisandro Dalcin /* Older versions of Chaco do not have an include file */
559371c9d4SSatish 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);
56abe9303eSLisandro Dalcin   #endif
57abe9303eSLisandro Dalcin extern int FREE_GRAPH;
58abe9303eSLisandro Dalcin #endif
59abe9303eSLisandro Dalcin 
6021c92275SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection edgeSection, PetscSection targetSection, PetscSection partSection, IS *partition)
61d71ae5a4SJacob Faibussowitsch {
62abe9303eSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
639371c9d4SSatish Balay   enum {
649371c9d4SSatish Balay     DEFAULT_METHOD  = 1,
659371c9d4SSatish Balay     INERTIAL_METHOD = 3
669371c9d4SSatish Balay   };
67abe9303eSLisandro Dalcin   MPI_Comm comm;
68abe9303eSLisandro Dalcin   int      nvtxs = numVertices;            /* number of vertices in full graph */
69abe9303eSLisandro Dalcin   int     *vwgts = NULL;                   /* weights for all vertices */
70abe9303eSLisandro Dalcin   float   *ewgts = NULL;                   /* weights for all edges */
71abe9303eSLisandro Dalcin   float   *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
72abe9303eSLisandro Dalcin   char    *outassignname = NULL;           /*  name of assignment output file */
73abe9303eSLisandro Dalcin   char    *outfilename   = NULL;           /* output file name */
74abe9303eSLisandro Dalcin   int      architecture  = 1;              /* 0 => hypercube, d => d-dimensional mesh */
75abe9303eSLisandro Dalcin   int      ndims_tot     = 0;              /* total number of cube dimensions to divide */
76abe9303eSLisandro Dalcin   int      mesh_dims[3];                   /* dimensions of mesh of processors */
77abe9303eSLisandro Dalcin   double  *goal          = NULL;           /* desired set sizes for each set */
78abe9303eSLisandro Dalcin   int      global_method = 1;              /* global partitioning algorithm */
79abe9303eSLisandro Dalcin   int      local_method  = 1;              /* local partitioning algorithm */
80abe9303eSLisandro Dalcin   int      rqi_flag      = 0;              /* should I use RQI/Symmlq eigensolver? */
81abe9303eSLisandro Dalcin   int      vmax          = 200;            /* how many vertices to coarsen down to? */
82abe9303eSLisandro Dalcin   int      ndims         = 1;              /* number of eigenvectors (2^d sets) */
83abe9303eSLisandro Dalcin   double   eigtol        = 0.001;          /* tolerance on eigenvectors */
84abe9303eSLisandro Dalcin   long     seed          = 123636512;      /* for random graph mutations */
85abe9303eSLisandro Dalcin   #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
86abe9303eSLisandro Dalcin   int *assignment; /* Output partition */
87abe9303eSLisandro Dalcin   #else
88abe9303eSLisandro Dalcin   short int *assignment; /* Output partition */
89abe9303eSLisandro Dalcin   #endif
90abe9303eSLisandro Dalcin   int       fd_stdout, fd_pipe[2];
91abe9303eSLisandro Dalcin   PetscInt *points;
92abe9303eSLisandro Dalcin   int       i, v, p;
93d0609cedSBarry Smith   int       err;
94abe9303eSLisandro Dalcin 
95abe9303eSLisandro Dalcin   PetscFunctionBegin;
969566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
97abe9303eSLisandro Dalcin   if (PetscDefined(USE_DEBUG)) {
98abe9303eSLisandro Dalcin     int       ival, isum;
99abe9303eSLisandro Dalcin     PetscBool distributed;
100abe9303eSLisandro Dalcin 
101abe9303eSLisandro Dalcin     ival = (numVertices > 0);
102*6a210b70SBarry Smith     PetscCallMPI(MPIU_Allreduce(&ival, &isum, 1, MPI_INT, MPI_SUM, comm));
103abe9303eSLisandro Dalcin     distributed = (isum > 1) ? PETSC_TRUE : PETSC_FALSE;
10428b400f6SJacob Faibussowitsch     PetscCheck(!distributed, comm, PETSC_ERR_SUP, "Chaco cannot partition a distributed graph");
105abe9303eSLisandro Dalcin   }
106abe9303eSLisandro Dalcin   if (!numVertices) { /* distributed case, return if not holding the graph */
1079566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition));
1083ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
109abe9303eSLisandro Dalcin   }
110abe9303eSLisandro Dalcin   FREE_GRAPH = 0; /* Do not let Chaco free my memory */
111abe9303eSLisandro Dalcin   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
112abe9303eSLisandro Dalcin 
113e08b1d6dSBarry Smith   /* code would use manager.createCellCoordinates(nvtxs, &x, &y, &z); */
114e08b1d6dSBarry Smith   PetscCheck(global_method != INERTIAL_METHOD, PETSC_COMM_SELF, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
115e08b1d6dSBarry Smith 
116abe9303eSLisandro Dalcin   mesh_dims[0] = nparts;
117abe9303eSLisandro Dalcin   mesh_dims[1] = 1;
118abe9303eSLisandro Dalcin   mesh_dims[2] = 1;
1199566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nvtxs, &assignment));
120abe9303eSLisandro Dalcin   /* Chaco outputs to stdout. We redirect this to a buffer. */
121abe9303eSLisandro Dalcin   /* TODO: check error codes for UNIX calls */
122abe9303eSLisandro Dalcin   #if defined(PETSC_HAVE_UNISTD_H)
123abe9303eSLisandro Dalcin   {
124abe9303eSLisandro Dalcin     int piperet;
125abe9303eSLisandro Dalcin     piperet = pipe(fd_pipe);
12628b400f6SJacob Faibussowitsch     PetscCheck(!piperet, PETSC_COMM_SELF, PETSC_ERR_SYS, "Could not create pipe");
127abe9303eSLisandro Dalcin     fd_stdout = dup(1);
128abe9303eSLisandro Dalcin     close(1);
129abe9303eSLisandro Dalcin     dup2(fd_pipe[1], 1);
130abe9303eSLisandro Dalcin   }
131abe9303eSLisandro Dalcin   #endif
1329566063dSJacob Faibussowitsch   if (part->usevwgt) PetscCall(PetscInfo(part, "PETSCPARTITIONERCHACO ignores vertex weights\n"));
13321c92275SMatthew G. Knepley   if (part->useewgt) PetscCall(PetscInfo(part, "PETSCPARTITIONERCHACO ignores edge weights\n"));
1349371c9d4SSatish 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);
135abe9303eSLisandro Dalcin   #if defined(PETSC_HAVE_UNISTD_H)
136abe9303eSLisandro Dalcin   {
137abe9303eSLisandro Dalcin     char msgLog[10000];
138abe9303eSLisandro Dalcin     int  count;
139abe9303eSLisandro Dalcin 
140c69effb2SJacob Faibussowitsch     PetscCall(PetscFFlush(stdout));
1416497c311SBarry Smith     count = (int)read(fd_pipe[0], msgLog, (10000 - 1) * sizeof(char));
142abe9303eSLisandro Dalcin     if (count < 0) count = 0;
143abe9303eSLisandro Dalcin     msgLog[count] = 0;
144abe9303eSLisandro Dalcin     close(1);
145abe9303eSLisandro Dalcin     dup2(fd_stdout, 1);
146abe9303eSLisandro Dalcin     close(fd_stdout);
147abe9303eSLisandro Dalcin     close(fd_pipe[0]);
148abe9303eSLisandro Dalcin     close(fd_pipe[1]);
149d0609cedSBarry Smith     PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
150abe9303eSLisandro Dalcin   }
151abe9303eSLisandro Dalcin   #else
152d0609cedSBarry Smith   PetscCheck(!err, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", "error in stdout");
153abe9303eSLisandro Dalcin   #endif
154abe9303eSLisandro Dalcin   /* Convert to PetscSection+IS */
15548a46eb9SPierre Jolivet   for (v = 0; v < nvtxs; ++v) PetscCall(PetscSectionAddDof(partSection, assignment[v], 1));
1569566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nvtxs, &points));
157abe9303eSLisandro Dalcin   for (p = 0, i = 0; p < nparts; ++p) {
158abe9303eSLisandro Dalcin     for (v = 0; v < nvtxs; ++v) {
159abe9303eSLisandro Dalcin       if (assignment[v] == p) points[i++] = v;
160abe9303eSLisandro Dalcin     }
161abe9303eSLisandro Dalcin   }
16263a3b9bcSJacob Faibussowitsch   PetscCheck(i == nvtxs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of points %" PetscInt_FMT " should be %" PetscInt_FMT, i, nvtxs);
1639566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition));
164e08b1d6dSBarry Smith   /* code would use manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
165e08b1d6dSBarry Smith   PetscCheck(global_method != INERTIAL_METHOD, PETSC_COMM_SELF, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
1669566063dSJacob Faibussowitsch   PetscCall(PetscFree(assignment));
167abe9303eSLisandro Dalcin   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
169abe9303eSLisandro Dalcin #else
170abe9303eSLisandro Dalcin   SETERRQ(PetscObjectComm((PetscObject)part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
171abe9303eSLisandro Dalcin #endif
172abe9303eSLisandro Dalcin }
173abe9303eSLisandro Dalcin 
174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
175d71ae5a4SJacob Faibussowitsch {
176abe9303eSLisandro Dalcin   PetscFunctionBegin;
177abe9303eSLisandro Dalcin   part->noGraph        = PETSC_FALSE;
178abe9303eSLisandro Dalcin   part->ops->view      = PetscPartitionerView_Chaco;
179abe9303eSLisandro Dalcin   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
180abe9303eSLisandro Dalcin   part->ops->partition = PetscPartitionerPartition_Chaco;
1813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182abe9303eSLisandro Dalcin }
183abe9303eSLisandro Dalcin 
184abe9303eSLisandro Dalcin /*MC
185abe9303eSLisandro Dalcin   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
186abe9303eSLisandro Dalcin 
187abe9303eSLisandro Dalcin   Level: intermediate
188abe9303eSLisandro Dalcin 
189db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
190abe9303eSLisandro Dalcin M*/
191abe9303eSLisandro Dalcin 
192d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
193d71ae5a4SJacob Faibussowitsch {
194abe9303eSLisandro Dalcin   PetscPartitioner_Chaco *p;
195abe9303eSLisandro Dalcin 
196abe9303eSLisandro Dalcin   PetscFunctionBegin;
197abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1984dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
199abe9303eSLisandro Dalcin   part->data = p;
200abe9303eSLisandro Dalcin 
2019566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Chaco(part));
2029566063dSJacob Faibussowitsch   PetscCall(PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionerCite));
2033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
204abe9303eSLisandro Dalcin }
205