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