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