Actual source code: chaco.c

 2:  #include src/mat/impls/adj/mpi/mpiadj.h

  4: #ifdef PETSC_HAVE_UNISTD_H
  5: #include <unistd.h>
  6: #endif

  8: #ifdef PETSC_HAVE_STDLIB_H
  9: #include <stdlib.h>
 10: #endif

 12: #include "petscfix.h"

 14: EXTERN_C_BEGIN
 15: /* Chaco does not have an include file */
 16: extern int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
 17:     float *ewgts, float *x, float *y, float *z, char *outassignname,
 18:     char *outfilename, short *assignment, int architecture, int ndims_tot,
 19:     int mesh_dims[3], double *goal, int global_method, int local_method,
 20:     int rqi_flag, int vmax, int ndims, double eigtol, long seed);

 22: extern int FREE_GRAPH;

 24: /*
 25: int       nvtxs;                number of vertices in full graph 
 26: int      *start;                start of edge list for each vertex 
 27: int      *adjacency;                edge list data 
 28: int      *vwgts;                weights for all vertices 
 29: float    *ewgts;                weights for all edges 
 30: float    *x, *y, *z;                coordinates for inertial method 
 31: char     *outassignname;        name of assignment output file 
 32: char     *outfilename;          output file name 
 33: short    *assignment;                set number of each vtx (length n) 
 34: int       architecture;         0 => hypercube, d => d-dimensional mesh 
 35: int       ndims_tot;                total number of cube dimensions to divide 
 36: int       mesh_dims[3];         dimensions of mesh of processors 
 37: double   *goal;                        desired set sizes for each set 
 38: int       global_method;        global partitioning algorithm 
 39: int       local_method;         local partitioning algorithm 
 40: int       rqi_flag;                should I use RQI/Symmlq eigensolver? 
 41: int       vmax;                        how many vertices to coarsen down to? 
 42: int       ndims;                number of eigenvectors (2^d sets) 
 43: double    eigtol;                tolerance on eigenvectors 
 44: long      seed;                        for random graph mutations 
 45: */

 47: EXTERN_C_END

 49: typedef struct {
 50:     int architecture;
 51:     int ndims_tot;
 52:     int mesh_dims[3];
 53:     int rqi_flag;
 54:     int numbereigen;
 55:     double eigtol;
 56:     int global_method;          /* global method */
 57:     int local_method;           /* local method */
 58:     int nbvtxcoarsed;           /* number of vertices for the coarse graph */
 59:     char *mesg_log;
 60: } MatPartitioning_Chaco;

 62: #define SIZE_LOG 10000          /* size of buffer for msg_log */

 66: static int MatPartitioningApply_Chaco(MatPartitioning part, IS *partitioning)
 67: {
 68:     int ierr, *parttab, *locals, i, size, rank;
 69:     Mat mat = part->adj, matMPI, matSeq;
 70:     int nb_locals;
 71:     Mat_MPIAdj *adj;
 72:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;
 73:     PetscTruth flg;
 74: #ifdef PETSC_HAVE_UNISTD_H
 75:     int fd_stdout, fd_pipe[2], count;
 76: #endif


 80:     FREE_GRAPH = 0; /* otherwise Chaco will attempt to free memory for adjacency graph */
 81: 
 82:     MPI_Comm_size(mat->comm, &size);

 84:     PetscTypeCompare((PetscObject) mat, MATMPIADJ, &flg);

 86:     /* check if the matrix is sequential, use MatGetSubMatrices if necessary */
 87:     if (size > 1) {
 88:         int M, N;
 89:         IS isrow, iscol;
 90:         Mat *A;

 92:         if (flg) {
 93:             SETERRQ(0, "Distributed matrix format MPIAdj is not supported for sequential partitioners");
 94:         }
 95:         PetscPrintf(part->comm, "Converting distributed matrix to sequential: this could be a performance loss\n");

 97:         MatGetSize(mat, &M, &N);
 98:         ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &isrow);
 99:         ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol);
100:         MatGetSubMatrices(mat, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &A);
101:         ISDestroy(isrow);
102:         ISDestroy(iscol);
103:         matSeq = *A;
104:     } else
105:         matSeq = mat;

107:     /* check for the input format that is supported only for a MPIADJ type 
108:        and set it to matMPI */
109:     if (!flg) {
110:         MatConvert(matSeq, MATMPIADJ, &matMPI);
111:     } else
112:         matMPI = matSeq;

114:     adj = (Mat_MPIAdj *) matMPI->data;  /* finaly adj contains adjacency graph */

116:     {
117:         /* arguments for Chaco library */
118:         int nvtxs = mat->M;                     /* number of vertices in full graph */
119:         int *start = adj->i;                    /* start of edge list for each vertex */
120:         int *adjacency;                         /* = adj -> j; edge list data  */
121:         int *vwgts = NULL;                      /* weights for all vertices */
122:         float *ewgts = NULL;                    /* weights for all edges */
123:         float *x = NULL, *y = NULL, *z = NULL;  /* coordinates for inertial method */
124:         char *outassignname = NULL;             /*  name of assignment output file */
125:         char *outfilename = NULL;               /* output file name */
126:         short *assignment;                      /* set number of each vtx (length n) */
127:         int architecture = chaco->architecture; /* 0 => hypercube, d => d-dimensional mesh */
128:         int ndims_tot = chaco->ndims_tot;       /* total number of cube dimensions to divide */
129:         int *mesh_dims = chaco->mesh_dims;      /* dimensions of mesh of processors */
130:         double *goal = NULL;                    /* desired set sizes for each set */
131:         int global_method = chaco->global_method; /* global partitioning algorithm */
132:         int local_method = chaco->local_method; /* local partitioning algorithm */
133:         int rqi_flag = chaco->rqi_flag;         /* should I use RQI/Symmlq eigensolver? */
134:         int vmax = chaco->nbvtxcoarsed;         /* how many vertices to coarsen down to? */
135:         int ndims = chaco->numbereigen;         /* number of eigenvectors (2^d sets) */
136:         double eigtol = chaco->eigtol;          /* tolerance on eigenvectors */
137:         long seed = 123636512;                  /* for random graph mutations */

139:         /* return value of Chaco */
140:         PetscMalloc((mat->M) * sizeof(short), &assignment);
141: 
142:         /* index change for libraries that have fortran implementation */
143:         PetscMalloc(sizeof(int) * start[nvtxs], &adjacency);
144:         for (i = 0; i < start[nvtxs]; i++)
145:             adjacency[i] = (adj->j)[i] + 1;

147:         /* redirect output to buffer: chaco -> mesg_log */
148: #ifdef PETSC_HAVE_UNISTD_H
149:         fd_stdout = dup(1);
150:         pipe(fd_pipe);
151:         close(1);
152:         dup2(fd_pipe[1], 1);
153:         PetscMalloc(SIZE_LOG * sizeof(char), &(chaco->mesg_log));
154: #endif

156:         /* library call */
157:         interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
158:             outassignname, outfilename, assignment, architecture, ndims_tot,
159:             mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
160:             eigtol, seed);

162: #ifdef PETSC_HAVE_UNISTD_H
163:         fflush(stdout);
164:         count =  read(fd_pipe[0], chaco->mesg_log, (SIZE_LOG - 1) * sizeof(char));
165:         if (count < 0)
166:             count = 0;
167:         chaco->mesg_log[count] = 0;
168:         close(1);
169:         dup2(fd_stdout, 1);
170:         close(fd_stdout);
171:         close(fd_pipe[0]);
172:         close(fd_pipe[1]);
173: #endif

175:         if (ierr) { SETERRQ(1, chaco->mesg_log); }

177:         PetscFree(adjacency);

179:         PetscMalloc((mat->M) * sizeof(int), &parttab);
180:         for (i = 0; i < nvtxs; i++) {
181:             parttab[i] = assignment[i];
182:         }
183:         PetscFree(assignment);
184:     }

186:     /* Creation of the index set */
187:     MPI_Comm_rank(part->comm, &rank);
188:     MPI_Comm_size(part->comm, &size);
189:     nb_locals = mat->M / size;
190:     locals = parttab + rank * nb_locals;
191:     if (rank < mat->M % size) {
192:         nb_locals++;
193:         locals += rank;
194:     } else
195:         locals += mat->M % size;
196:     ISCreateGeneral(part->comm, nb_locals, locals, partitioning);

198:     /* destroy temporary objects */
199:     PetscFree(parttab);
200:     if (matSeq != mat) {
201:         MatDestroy(matSeq);
202:     }
203:     if (matMPI != mat) {
204:         MatDestroy(matMPI);
205:     }

207:     return(0);
208: }


213: int MatPartitioningView_Chaco(MatPartitioning part, PetscViewer viewer)
214: {

216:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;
217:     int ierr, rank;
218:     PetscTruth isascii;


222:     MPI_Comm_rank(part->comm, &rank);
223:     PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &isascii);
224:     if (isascii) {
225:         if (!rank && chaco->mesg_log) {
226:             PetscViewerASCIIPrintf(viewer, "%s\n", chaco->mesg_log);
227:         }
228:     } else {
229:         SETERRQ1(1, "Viewer type %s not supported for this Chaco partitioner",
230:             ((PetscObject) viewer)->type_name);
231:     }

233:     return(0);
234: }

238: /*@
239:      MatPartitioningChacoSetGlobal - Set method for global partitioning.

241:   Input Parameter:
242: .  part - the partitioning context
243: .  method - MP_CHACO_MULTILEVEL_KL, MP_CHACO_SPECTRAL, MP_CHACO_LINEAR, 
244:     MP_CHACO_RANDOM or MP_CHACO_SCATTERED

246:    Level: advanced

248: @*/
249: int MatPartitioningChacoSetGlobal(MatPartitioning part, MPChacoGlobalType method)
250: {
251:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


255:     switch (method) {
256:     case MP_CHACO_MULTILEVEL_KL:
257:         chaco->global_method = 1;
258:         break;
259:     case MP_CHACO_SPECTRAL:
260:         chaco->global_method = 2;
261:         break;
262:     case MP_CHACO_LINEAR:
263:         chaco->global_method = 4;
264:         break;
265:     case MP_CHACO_RANDOM:
266:         chaco->global_method = 5;
267:         break;
268:     case MP_CHACO_SCATTERED:
269:         chaco->global_method = 6;
270:         break;
271:     default:
272:         SETERRQ(PETSC_ERR_SUP, "Chaco: Unknown or unsupported option");
273:     }
274:     return(0);
275: }

279: /*@
280:      MatPartitioningChacoSetLocal - Set method for local partitioning.

282:   Input Parameter:
283: .  part - the partitioning context
284: .  method - MP_CHACO_KERNIGHAN_LIN or MP_CHACO_NONE

286:    Level: advanced

288: @*/
289: int MatPartitioningChacoSetLocal(MatPartitioning part, MPChacoLocalType method)
290: {
291:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


295:     switch (method) {
296:     case MP_CHACO_KERNIGHAN_LIN:
297:         chaco->local_method = 1;
298:         break;
299:     case MP_CHACO_NONE:
300:         chaco->local_method = 2;
301:         break;
302:     default:
303:         SETERRQ(PETSC_ERR_ARG_CORRUPT, "Chaco: Unknown or unsupported option");
304:     }

306:     return(0);
307: }

311: /*@
312:     MatPartitioningChacoSetCoarseLevel - Set the coarse level 
313:     
314:   Input Parameter:
315: .  part - the partitioning context
316: .  level - the coarse level in range [0.0,1.0]

318:    Level: advanced

320: @*/
321: int MatPartitioningChacoSetCoarseLevel(MatPartitioning part, PetscReal level)
322: {
323:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


327:     if (level < 0 || level > 1.0) {
328:         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,
329:             "Chaco: level of coarsening out of range [0.01-1.0]");
330:     } else
331:         chaco->nbvtxcoarsed = part->adj->N * level;

333:     if (chaco->nbvtxcoarsed < 20)
334:         chaco->nbvtxcoarsed = 20;

336:     return(0);
337: }

341: /*@
342:      MatPartitioningChacoSetEigenSolver - Set method for eigensolver.

344:   Input Parameter:
345: .  method - MP_CHACO_LANCZOS or MP_CHACO_RQI_SYMMLQ

347:    Level: advanced

349: @*/
350: int MatPartitioningChacoSetEigenSolver(MatPartitioning part,
351:     MPChacoEigenType method)
352: {
353:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


357:     switch (method) {
358:     case MP_CHACO_LANCZOS:
359:         chaco->rqi_flag = 0;
360:         break;
361:     case MP_CHACO_RQI_SYMMLQ:
362:         chaco->rqi_flag = 1;
363:         break;
364:     default:
365:         SETERRQ(PETSC_ERR_ARG_CORRUPT, "Chaco: Unknown or unsupported option");
366:     }

368:     return(0);
369: }

373: /*@
374:      MatPartitioningChacoSetEigenTol - Set tolerance for eigensolver.

376:   Input Parameter:
377: .  tol - Tolerance requested.

379:    Level: advanced

381: @*/
382: int MatPartitioningChacoSetEigenTol(MatPartitioning part, PetscReal tol)
383: {
384:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


388:     if (tol <= 0.0) {
389:         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,
390:             "Chaco: Eigensolver tolerance out of range");
391:     } else
392:         chaco->eigtol = tol;

394:     return(0);
395: }

399: /*@
400:      MatPartitioningChacoSetEigenNumber - Set number of eigenvectors for partitioning.

402:   Input Parameter:
403: .  num - This argument should have a value of 1, 2 or 3 indicating  
404:     partitioning by bisection, quadrisection, or octosection.

406:    Level: advanced

408: @*/
409: int MatPartitioningChacoSetEigenNumber(MatPartitioning part, int num)
410: {
411:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


415:     if (num > 3 || num < 1) {
416:         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,
417:             "Chaco: number of eigenvectors out of range");
418:     } else
419:         chaco->numbereigen = num;

421:     return(0);
422: }

426: int MatPartitioningSetFromOptions_Chaco(MatPartitioning part)
427: {
428:     int ierr, i;
429:     PetscReal r;
430:     PetscTruth flag;
431:     const char *global[] =
432:         { "multilevel-kl", "spectral", "linear", "random", "scattered" };
433:     const char *local[] = { "kernighan-lin", "none" };
434:     const char *eigen[] = { "lanczos", "rqi_symmlq" };

437:     PetscOptionsHead("Set Chaco partitioning options");

439:     PetscOptionsEList("-mat_partitioning_chaco_global",
440:         "Global method to use", "MatPartitioningChacoSetGlobal", global, 5,
441:         global[0], &i, &flag);
442:     if (flag)
443:         MatPartitioningChacoSetGlobal(part, i);

445:     PetscOptionsEList("-mat_partitioning_chaco_local",
446:         "Local method to use", "MatPartitioningChacoSetLocal", local, 2,
447:         local[0], &i, &flag);
448:     if (flag)
449:         MatPartitioningChacoSetLocal(part, i);

451:     PetscOptionsReal("-mat_partitioning_chaco_coarse_level",
452:         "Coarse level", "MatPartitioningChacoSetCoarseLevel", 0, &r,
453:         &flag);
454:     if (flag)
455:         MatPartitioningChacoSetCoarseLevel(part, r);

457:     PetscOptionsEList("-mat_partitioning_chaco_eigen_solver",
458:         "Eigensolver to use in spectral method", "MatPartitioningChacoSetEigenSolver",
459:         eigen, 2, eigen[0], &i, &flag);
460:     if (flag)
461:         MatPartitioningChacoSetEigenSolver(part, i);

463:     PetscOptionsReal("-mat_partitioning_chaco_eigen_tol",
464:         "Tolerance for eigensolver", "MatPartitioningChacoSetEigenTol", 0.001,
465:         &r, &flag);
466:     if (flag)
467:         MatPartitioningChacoSetEigenTol(part, r);

469:     PetscOptionsInt("-mat_partitioning_chaco_eigen_number",
470:         "Number of eigenvectors: 1, 2, or 3 (bi-, quadri-, or octosection)",
471:         "MatPartitioningChacoSetEigenNumber", 1, &i, &flag);
472:     if (flag)
473:         MatPartitioningChacoSetEigenNumber(part, i);

475:     PetscOptionsTail();
476:     return(0);
477: }


482: int MatPartitioningDestroy_Chaco(MatPartitioning part)
483: {
484:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


489:     if (chaco->mesg_log) {
490:         PetscFree(chaco->mesg_log);
491:     }

493:     PetscFree(chaco);

495:     return(0);
496: }

498: EXTERN_C_BEGIN
501: int MatPartitioningCreate_Chaco(MatPartitioning part)
502: {
504:     MatPartitioning_Chaco *chaco;

507:     PetscNew(MatPartitioning_Chaco, &chaco);

509:     chaco->architecture = 1;
510:     chaco->ndims_tot = 0;
511:     chaco->mesh_dims[0] = part->n;
512:     chaco->global_method = 1;
513:     chaco->local_method = 1;
514:     chaco->rqi_flag = 0;
515:     chaco->nbvtxcoarsed = 200;
516:     chaco->numbereigen = 1;
517:     chaco->eigtol = 0.001;
518:     chaco->mesg_log = NULL;

520:     part->ops->apply = MatPartitioningApply_Chaco;
521:     part->ops->view = MatPartitioningView_Chaco;
522:     part->ops->destroy = MatPartitioningDestroy_Chaco;
523:     part->ops->setfromoptions = MatPartitioningSetFromOptions_Chaco;
524:     part->data = (void *) chaco;

526:     return(0);
527: }

529: EXTERN_C_END