Actual source code: chaco.c

  1: #define PETSCMAT_DLL

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

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

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

 13: #include "petscfix.h"

 16: /* Chaco does not have an include file */
 18:     float *ewgts, float *x, float *y, float *z, char *outassignname,
 19:     char *outfilename, short *assignment, int architecture, int ndims_tot,
 20:     int mesh_dims[3], double *goal, int global_method, int local_method,
 21:     int rqi_flag, int vmax, int ndims, double eigtol, long seed);


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


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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

178:         if (ierr) { SETERRQ(PETSC_ERR_LIB, chaco->mesg_log); }

180:         PetscFree(adjacency);

182:         PetscMalloc((mat->rmap.N) * sizeof(int), &parttab);
183:         for (i = 0; i < nvtxs; i++) {
184:             parttab[i] = assignment[i];
185:         }
186:         PetscFree(assignment);
187:     }

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

201:     /* destroy temporary objects */
202:     PetscFree(parttab);
203:     if (matSeq != mat) {
204:         MatDestroy(matSeq);
205:     }
206:     if (matMPI != mat) {
207:         MatDestroy(matMPI);
208:     }
209:     return(0);
210: }


215: PetscErrorCode MatPartitioningView_Chaco(MatPartitioning part, PetscViewer viewer)
216: {

218:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;
219:     PetscErrorCode        ierr;
220:     PetscMPIInt           rank;
221:     PetscTruth            iascii;


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

235:     return(0);
236: }

240: /*@
241:      MatPartitioningChacoSetGlobal - Set method for global partitioning.

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

248:    Level: advanced

250: @*/
251: PetscErrorCode  MatPartitioningChacoSetGlobal(MatPartitioning part, MPChacoGlobalType method)
252: {
253:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


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

281: /*@
282:      MatPartitioningChacoSetLocal - Set method for local partitioning.

284:   Input Parameter:
285: .  part - the partitioning context
286: .  method - MP_CHACO_KERNIGHAN_LIN or MP_CHACO_NONE

288:    Level: advanced

290: @*/
291: PetscErrorCode  MatPartitioningChacoSetLocal(MatPartitioning part, MPChacoLocalType method)
292: {
293:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


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

308:     return(0);
309: }

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

320:    Level: advanced

322: @*/
323: PetscErrorCode  MatPartitioningChacoSetCoarseLevel(MatPartitioning part, PetscReal level)
324: {
325:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


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

335:     if (chaco->nbvtxcoarsed < 20)
336:         chaco->nbvtxcoarsed = 20;

338:     return(0);
339: }

343: /*@
344:      MatPartitioningChacoSetEigenSolver - Set method for eigensolver.

346:   Input Parameter:
347: .  method - MP_CHACO_LANCZOS or MP_CHACO_RQI_SYMMLQ

349:    Level: advanced

351: @*/
352: PetscErrorCode  MatPartitioningChacoSetEigenSolver(MatPartitioning part,
353:     MPChacoEigenType method)
354: {
355:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


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

370:     return(0);
371: }

375: /*@
376:      MatPartitioningChacoSetEigenTol - Set tolerance for eigensolver.

378:   Input Parameter:
379: .  tol - Tolerance requested.

381:    Level: advanced

383: @*/
384: PetscErrorCode  MatPartitioningChacoSetEigenTol(MatPartitioning part, PetscReal tol)
385: {
386:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;


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

396:     return(0);
397: }

401: /*@
402:      MatPartitioningChacoSetEigenNumber - Set number of eigenvectors for partitioning.

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

408:    Level: advanced

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


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

423:     return(0);
424: }

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

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

442:     PetscOptionsEList("-mat_partitioning_chaco_global",
443:         "Global method to use", "MatPartitioningChacoSetGlobal", global, 5,
444:         global[0], &i, &flag);
445:     if (flag)
446:         MatPartitioningChacoSetGlobal(part, (MPChacoGlobalType)i);

448:     PetscOptionsEList("-mat_partitioning_chaco_local",
449:         "Local method to use", "MatPartitioningChacoSetLocal", local, 2,
450:         local[0], &i, &flag);
451:     if (flag)
452:         MatPartitioningChacoSetLocal(part, (MPChacoLocalType)i);

454:     PetscOptionsReal("-mat_partitioning_chaco_coarse_level",
455:         "Coarse level", "MatPartitioningChacoSetCoarseLevel", 0, &r,
456:         &flag);
457:     if (flag)
458:         MatPartitioningChacoSetCoarseLevel(part, r);

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

466:     PetscOptionsReal("-mat_partitioning_chaco_eigen_tol",
467:         "Tolerance for eigensolver", "MatPartitioningChacoSetEigenTol", 0.001,
468:         &r, &flag);
469:     if (flag)
470:         MatPartitioningChacoSetEigenTol(part, r);

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

478:     PetscOptionsTail();
479:     return(0);
480: }


485: PetscErrorCode MatPartitioningDestroy_Chaco(MatPartitioning part)
486: {
487:     MatPartitioning_Chaco *chaco = (MatPartitioning_Chaco *) part->data;
488:     PetscErrorCode        ierr;

491:     PetscFree(chaco->mesg_log);
492:     PetscFree(chaco);
493:     return(0);
494: }

499: /*@C
500:    MAT_PARTITIONING_Chaco - Creates a partitioning context via the external package Chaco.

502:    Collective on MPI_Comm

504:    Input Parameter:
505: .  part - the partitioning context

507:    Options Database Keys:
508: +  -mat_partitioning_chaco_global <multilevel-kl> (one of) multilevel-kl spectral linear random scattered
509: .  -mat_partitioning_chaco_local <kernighan-lin> (one of) kernighan-lin none
510: .  -mat_partitioning_chaco_coarse_level <0>: Coarse level (MatPartitioningChacoSetCoarseLevel)
511: .  -mat_partitioning_chaco_eigen_solver <lanczos> (one of) lanczos rqi_symmlq
512: .  -mat_partitioning_chaco_eigen_tol <0.001>: Tolerance for eigensolver (MatPartitioningChacoSetEigenTol)
513: -  -mat_partitioning_chaco_eigen_number <1>: Number of eigenvectors: 1, 2, or 3 (bi-, quadri-, or octosection) (MatPartitioningChacoSetEigenNumber)

515:    Level: beginner

517:    Notes: See http://www.cs.sandia.gov/CRF/chac.html

519: .keywords: Partitioning, create, context

521: .seealso: MatPartitioningSetType(), MatPartitioningType

523: @*/
524: PetscErrorCode  MatPartitioningCreate_Chaco(MatPartitioning part)
525: {
527:     MatPartitioning_Chaco *chaco;

530:     PetscNew(MatPartitioning_Chaco, &chaco);

532:     chaco->architecture = 1;
533:     chaco->ndims_tot = 0;
534:     chaco->mesh_dims[0] = part->n;
535:     chaco->global_method = 1;
536:     chaco->local_method = 1;
537:     chaco->rqi_flag = 0;
538:     chaco->nbvtxcoarsed = 200;
539:     chaco->numbereigen = 1;
540:     chaco->eigtol = 0.001;
541:     chaco->mesg_log = NULL;

543:     part->ops->apply = MatPartitioningApply_Chaco;
544:     part->ops->view = MatPartitioningView_Chaco;
545:     part->ops->destroy = MatPartitioningDestroy_Chaco;
546:     part->ops->setfromoptions = MatPartitioningSetFromOptions_Chaco;
547:     part->data = (void*) chaco;

549:     return(0);
550: }