Actual source code: partchaco.c

  1: #include <petsc/private/partitionerimpl.h>

  3: PetscBool  ChacoPartitionerCite = PETSC_FALSE;
  4: const char ChacoPartitionerCitation[] =
  5:   "@inproceedings{Chaco95,\n"
  6:   "  author    = {Bruce Hendrickson and Robert Leland},\n"
  7:   "  title     = {A multilevel algorithm for partitioning graphs},\n"
  8:   "  booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)},"
  9:   "  isbn      = {0-89791-816-9},\n"
 10:   "  pages     = {28},\n"
 11:   "  doi       = {https://doi.acm.org/10.1145/224170.224228},\n"
 12:   "  publisher = {ACM Press},\n"
 13:   "  address   = {New York},\n"
 14:   "  year      = {1995}\n"
 15:   "}\n";

 17: typedef struct {
 18:   PetscInt dummy;
 19: } PetscPartitioner_Chaco;

 21: static PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
 22: {
 23:   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data;
 24:   PetscErrorCode          ierr;

 27:   PetscFree(p);
 28:   return(0);
 29: }

 31: static PetscErrorCode PetscPartitionerView_Chaco_ASCII(PetscPartitioner part, PetscViewer viewer)
 32: {
 34:   return(0);
 35: }

 37: static PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
 38: {
 39:   PetscBool      iascii;

 45:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);
 46:   if (iascii) {PetscPartitionerView_Chaco_ASCII(part, viewer);}
 47:   return(0);
 48: }

 50: #if defined(PETSC_HAVE_CHACO)
 51: #if defined(PETSC_HAVE_UNISTD_H)
 52: #include <unistd.h>
 53: #endif
 54: #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
 55: #include <chaco.h>
 56: #else
 57: /* Older versions of Chaco do not have an include file */
 58: PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
 59:                        float *ewgts, float *x, float *y, float *z, char *outassignname,
 60:                        char *outfilename, short *assignment, int architecture, int ndims_tot,
 61:                        int mesh_dims[3], double *goal, int global_method, int local_method,
 62:                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
 63: #endif
 64: extern int FREE_GRAPH;
 65: #endif

 67: static PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
 68: {
 69: #if defined(PETSC_HAVE_CHACO)
 70:   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
 71:   MPI_Comm       comm;
 72:   int            nvtxs          = numVertices; /* number of vertices in full graph */
 73:   int           *vwgts          = NULL;   /* weights for all vertices */
 74:   float         *ewgts          = NULL;   /* weights for all edges */
 75:   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
 76:   char          *outassignname  = NULL;   /*  name of assignment output file */
 77:   char          *outfilename    = NULL;   /* output file name */
 78:   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
 79:   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
 80:   int            mesh_dims[3];            /* dimensions of mesh of processors */
 81:   double        *goal          = NULL;    /* desired set sizes for each set */
 82:   int            global_method = 1;       /* global partitioning algorithm */
 83:   int            local_method  = 1;       /* local partitioning algorithm */
 84:   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
 85:   int            vmax          = 200;     /* how many vertices to coarsen down to? */
 86:   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
 87:   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
 88:   long           seed          = 123636512; /* for random graph mutations */
 89: #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
 90:   int           *assignment;              /* Output partition */
 91: #else
 92:   short int     *assignment;              /* Output partition */
 93: #endif
 94:   int            fd_stdout, fd_pipe[2];
 95:   PetscInt      *points;
 96:   int            i, v, p;

100:   PetscObjectGetComm((PetscObject)part,&comm);
101:   if (PetscDefined (USE_DEBUG)) {
102:     int ival,isum;
103:     PetscBool distributed;

105:     ival = (numVertices > 0);
106:     MPI_Allreduce(&ival, &isum, 1, MPI_INT, MPI_SUM, comm);
107:     distributed = (isum > 1) ? PETSC_TRUE : PETSC_FALSE;
108:     if (distributed) SETERRQ(comm, PETSC_ERR_SUP, "Chaco cannot partition a distributed graph");
109:   }
110:   if (!numVertices) { /* distributed case, return if not holding the graph */
111:     ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);
112:     return(0);
113:   }
114:   FREE_GRAPH = 0; /* Do not let Chaco free my memory */
115:   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];

117:   if (global_method == INERTIAL_METHOD) {
118:     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
119:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
120:   }
121:   mesh_dims[0] = nparts;
122:   mesh_dims[1] = 1;
123:   mesh_dims[2] = 1;
124:   PetscMalloc1(nvtxs, &assignment);
125:   /* Chaco outputs to stdout. We redirect this to a buffer. */
126:   /* TODO: check error codes for UNIX calls */
127: #if defined(PETSC_HAVE_UNISTD_H)
128:   {
129:     int piperet;
130:     piperet = pipe(fd_pipe);
131:     if (piperet) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Could not create pipe");
132:     fd_stdout = dup(1);
133:     close(1);
134:     dup2(fd_pipe[1], 1);
135:   }
136: #endif
137:   if (part->usevwgt) { PetscInfo(part,"PETSCPARTITIONERCHACO ignores vertex weights\n"); }
138:   interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
139:                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
140:                    vmax, ndims, eigtol, seed);
141: #if defined(PETSC_HAVE_UNISTD_H)
142:   {
143:     char msgLog[10000];
144:     int  count;

146:     fflush(stdout);
147:     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
148:     if (count < 0) count = 0;
149:     msgLog[count] = 0;
150:     close(1);
151:     dup2(fd_stdout, 1);
152:     close(fd_stdout);
153:     close(fd_pipe[0]);
154:     close(fd_pipe[1]);
155:     if (ierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
156:   }
157: #else
158:   if (ierr) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in Chaco library: %s", "error in stdout");
159: #endif
160:   /* Convert to PetscSection+IS */
161:   for (v = 0; v < nvtxs; ++v) {
162:     PetscSectionAddDof(partSection, assignment[v], 1);
163:   }
164:   PetscMalloc1(nvtxs, &points);
165:   for (p = 0, i = 0; p < nparts; ++p) {
166:     for (v = 0; v < nvtxs; ++v) {
167:       if (assignment[v] == p) points[i++] = v;
168:     }
169:   }
170:   if (i != nvtxs) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
171:   ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);
172:   if (global_method == INERTIAL_METHOD) {
173:     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
174:   }
175:   PetscFree(assignment);
176:   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
177:   return(0);
178: #else
179:   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
180: #endif
181: }

183: static PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
184: {
186:   part->noGraph        = PETSC_FALSE;
187:   part->ops->view      = PetscPartitionerView_Chaco;
188:   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
189:   part->ops->partition = PetscPartitionerPartition_Chaco;
190:   return(0);
191: }

193: /*MC
194:   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library

196:   Level: intermediate

198: .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
199: M*/

201: PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
202: {
203:   PetscPartitioner_Chaco *p;
204:   PetscErrorCode          ierr;

208:   PetscNewLog(part, &p);
209:   part->data = p;

211:   PetscPartitionerInitialize_Chaco(part);
212:   PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionerCite);
213:   return(0);
214: }