Actual source code: scotch.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"

 15: /* 
 16:    Currently using Scotch-3.4
 17: */
 19: #include "scotch.h"

 22: typedef struct {
 23:     char arch[PETSC_MAX_PATH_LEN];
 24:     int multilevel;
 25:     char strategy[30];
 26:     int global_method;          /* global method */
 27:     int local_method;           /* local method */
 28:     int nbvtxcoarsed;           /* number of vertices for the coarse graph */
 29:     int map;                    /* to know if we map on archptr or just partionate the graph */
 30:     char *mesg_log;
 31:     char host_list[PETSC_MAX_PATH_LEN];
 32: } MatPartitioning_Scotch;

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

 38: static PetscErrorCode MatPartitioningApply_Scotch(MatPartitioning part, IS * partitioning)
 39: {
 41:     int  *parttab, *locals = PETSC_NULL, rank, i, size;
 42:     size_t                 j;
 43:     Mat                    mat = part->adj, matMPI, matSeq;
 44:     int                    nb_locals = mat->m;
 45:     Mat_MPIAdj             *adj = (Mat_MPIAdj *) mat->data;
 46:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;
 47:     PetscTruth             flg;
 48: #ifdef PETSC_HAVE_UNISTD_H
 49:     int                    fd_stdout, fd_pipe[2], count;
 50: #endif


 54:     /* check if the matrix is sequential, use MatGetSubMatrices if necessary */
 55:     MPI_Comm_size(mat->comm, &size);
 56:     PetscTypeCompare((PetscObject) mat, MATMPIADJ, &flg);
 57:     if (size > 1) {
 58:         int M, N;
 59:         IS isrow, iscol;
 60:         Mat *A;

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

 67:         MatGetSize(mat, &M, &N);
 68:         ISCreateStride(PETSC_COMM_SELF, M, 0, 1, &isrow);
 69:         ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &iscol);
 70:         MatGetSubMatrices(mat, 1, &isrow, &iscol, MAT_INITIAL_MATRIX, &A);
 71:         matSeq = *A;
 72:         PetscFree(A);
 73:         ISDestroy(isrow);
 74:         ISDestroy(iscol);
 75:     } else
 76:         matSeq = mat;

 78:     /* convert the the matrix to MPIADJ type if necessary */
 79:     if (!flg) {
 80:         MatConvert(matSeq, MATMPIADJ, MAT_INITIAL_MATRIX, &matMPI);
 81:     } else
 82:         matMPI = matSeq;

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

 86:     MPI_Comm_rank(part->comm, &rank);

 88:     {
 89:         /* definition of Scotch library arguments */
 90:         SCOTCH_Strat stratptr;  /* scotch strategy */
 91:         SCOTCH_Graph grafptr;   /* scotch graph */
 92:         SCOTCH_Mapping mappptr; /* scotch mapping format */
 93:         int vertnbr = mat->M;   /* number of vertices in full graph */
 94:         int *verttab = adj->i;  /* start of edge list for each vertex */
 95:         int *edgetab = adj->j;  /* edge list data */
 96:         int edgenbr = adj->nz;  /* number of edges */
 97:         int *velotab = NULL;    /* not used by petsc interface */
 98:         int *vlbltab = NULL;
 99:         int *edlotab = NULL;
100:         int baseval = 0;        /* 0 for C array indexing */
101:         int flagval = 3;        /* (cf doc scotch no weight edge & vertices) */
102:         char strategy[256];

104:         PetscMalloc((mat->M) * sizeof(int), &parttab);

106:         /* redirect output to buffer scotch -> mesg_log */
107: #ifdef PETSC_HAVE_UNISTD_H
108:         fd_stdout = dup(1);
109:         pipe(fd_pipe);
110:         close(1);
111:         dup2(fd_pipe[1], 1);
112:         PetscMalloc(SIZE_LOG * sizeof(char), &(scotch->mesg_log));
113: #endif

115:         /* library call */

117:         /* Construction of the scotch graph object */
118:         SCOTCH_graphInit(&grafptr);
119:         SCOTCH_graphBuild(&grafptr, vertnbr, verttab, velotab,
120:             vlbltab, edgenbr, edgetab, edlotab, baseval, flagval);
121:         SCOTCH_graphCheck(&grafptr);

123:         /* Construction of the strategy */
124:         if (scotch->strategy[0] != 0)   /* strcmp(scotch->strategy,"") */
125:             PetscStrcpy(strategy, scotch->strategy);
126:         else {
127:             PetscStrcpy(strategy, "b{strat=");

129:             if (scotch->multilevel) {
130:                 /* PetscStrcat(strategy,"m{vert=");
131:                    sprintf(strategy+strlen(strategy),"%d",scotch->nbvtxcoarsed);
132:                    PetscStrcat(strategy,",asc="); */
133:                 sprintf(strategy, "b{strat=m{vert=%d,asc=",
134:                     scotch->nbvtxcoarsed);
135:             } else
136:                 PetscStrcpy(strategy, "b{strat=");

138:             switch (scotch->global_method) {
139:             case MP_SCOTCH_GREEDY:
140:                 PetscStrcat(strategy, "h");
141:                 break;
142:             case MP_SCOTCH_GPS:
143:                 PetscStrcat(strategy, "g");
144:                 break;
145:             case MP_SCOTCH_GR_GPS:
146:                 PetscStrcat(strategy, "g|h");
147:             }

149:             switch (scotch->local_method) {
150:             case MP_SCOTCH_KERNIGHAN_LIN:
151:                 if (scotch->multilevel)
152:                     PetscStrcat(strategy, ",low=f}");
153:                 else
154:                     PetscStrcat(strategy, " f");
155:                 break;
156:             case MP_SCOTCH_NONE:
157:                 if (scotch->multilevel)
158:                     PetscStrcat(strategy, ",asc=x}");
159:             default:
160:                 break;
161:             }

163:             PetscStrcat(strategy, " x}");
164:         }

166:         PetscPrintf(part->comm, "strategy=[%s]\n", strategy);

168:         SCOTCH_stratInit(&stratptr);
169:         SCOTCH_stratMap(&stratptr, strategy);

171:         /* check for option mapping */
172:         if (!scotch->map) {
173:             SCOTCH_graphPart(&grafptr, &stratptr, part->n, parttab);
174:             PetscPrintf(PETSC_COMM_SELF, "Partition simple without mapping\n");
175:         } else {
176:             SCOTCH_Graph grafarch;
177:             SCOTCH_Num *listtab;
178:             SCOTCH_Num listnbr = 0;
179:             SCOTCH_Arch archptr;        /* file in scotch architecture format */
180:             SCOTCH_Strat archstrat;
181:             int arch_total_size, *parttab_tmp;
182:             int cpt;
183:             char buf[256];
184:             FILE *file1, *file2;
185:             char host_buf[256];

187:             /* generate the graph that represents the arch */
188:             file1 = fopen(scotch->arch, "r");
189:             if (!file1)
190:                 SETERRQ1(PETSC_ERR_FILE_OPEN, "Scotch: unable to open architecture file %s", scotch->arch);

192:             SCOTCH_graphInit(&grafarch);
193:             SCOTCH_graphLoad(&grafarch, file1, baseval, 3);

195:             SCOTCH_graphCheck(&grafarch);
196:             SCOTCH_graphSize(&grafarch, &arch_total_size, &cpt);

198:             fclose(file1);
199:             printf("total size = %d\n", arch_total_size);

201:             /* generate the list of nodes currently working */
202:             PetscGetHostName(host_buf, 256);
203:             PetscStrlen(host_buf, &j);

205:             file2 = fopen(scotch->host_list, "r");
206:             if (!file2)
207:                 SETERRQ1(PETSC_ERR_FILE_OPEN, "Scotch: unable to open host list file %s", scotch->host_list);

209:             i = -1;
210:             flg = PETSC_FALSE;
211:             while (!feof(file2) && !flg) {
212:                 i++;
213:                 fgets(buf, 256, file2);
214:                 PetscStrncmp(buf, host_buf, j, &flg);
215:             }
216:             fclose(file2);
217:             if (!flg) {
218:                 SETERRQ1(PETSC_ERR_LIB, "Scotch: unable to find '%s' in host list file", host_buf);
219:             }

221:             listnbr = size;
222:             PetscMalloc(sizeof(SCOTCH_Num) * listnbr, &listtab);

224:             MPI_Allgather(&i, 1, MPI_INT, listtab, 1, MPI_INT, part->comm);

226:             printf("listnbr = %d, listtab = ", listnbr);
227:             for (i = 0; i < listnbr; i++)
228:                 printf("%d ", listtab[i]);

230:             printf("\n");
231:             fflush(stdout);

233:             SCOTCH_stratInit(&archstrat);
234:             SCOTCH_stratBipart(&archstrat, "fx");

236:             SCOTCH_archInit(&archptr);
237:             SCOTCH_archBuild(&archptr, &grafarch, listnbr, listtab,
238:                 &archstrat);

240:             PetscMalloc((mat->M) * sizeof(int), &parttab_tmp);
241:             SCOTCH_mapInit(&mappptr, &grafptr, &archptr, parttab_tmp);

243:             SCOTCH_mapCompute(&mappptr, &stratptr);

245:             SCOTCH_mapView(&mappptr, stdout);

247:             /* now we have to set in the real parttab at the good place */
248:             /* because the ranks order are different than position in */
249:             /* the arch graph */
250:             for (i = 0; i < mat->M; i++) {
251:                 parttab[i] = parttab_tmp[i];
252:             }

254:             PetscFree(listtab);
255:             SCOTCH_archExit(&archptr);
256:             SCOTCH_mapExit(&mappptr);
257:             SCOTCH_stratExit(&archstrat);
258:         }

260:         /* dump to mesg_log... */
261: #ifdef PETSC_HAVE_UNISTD_H
262:         fflush(stdout);
263:         count = read(fd_pipe[0], scotch->mesg_log, (SIZE_LOG - 1) * sizeof(char));
264:         if (count < 0)
265:             count = 0;
266:         scotch->mesg_log[count] = 0;
267:         close(1);
268:         dup2(fd_stdout, 1);
269:         close(fd_stdout);
270:         close(fd_pipe[0]);
271:         close(fd_pipe[1]);
272: #endif

274:         SCOTCH_graphExit(&grafptr);
275:         SCOTCH_stratExit(&stratptr);
276:     }

278:     if (ierr)
279:         SETERRQ(PETSC_ERR_LIB, scotch->mesg_log);

281:     /* Creation of the index set */

283:     MPI_Comm_rank(part->comm, &rank);
284:     MPI_Comm_size(part->comm, &size);
285:     nb_locals = mat->M / size;
286:     locals = parttab + rank * nb_locals;
287:     if (rank < mat->M % size) {
288:         nb_locals++;
289:         locals += rank;
290:     } else
291:         locals += mat->M % size;
292:     ISCreateGeneral(part->comm, nb_locals, locals, partitioning);

294:     /* destroying old objects */
295:     PetscFree(parttab);
296:     if (matSeq != mat) {
297:         MatDestroy(matSeq);
298:     }
299:     if (matMPI != mat) {
300:         MatDestroy(matMPI);
301:     }

303:     return(0);
304: }


309: PetscErrorCode MatPartitioningView_Scotch(MatPartitioning part, PetscViewer viewer)
310: {
311:   MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;
312:   PetscErrorCode         ierr;
313:   PetscMPIInt            rank;
314:   PetscTruth             iascii;
315: 
317:   MPI_Comm_rank(part->comm, &rank);
318:   PetscTypeCompare((PetscObject) viewer, PETSC_VIEWER_ASCII, &iascii);
319:   if (iascii) {
320:     if (!rank && scotch->mesg_log) {
321:       PetscViewerASCIIPrintf(viewer, "%s\n", scotch->mesg_log);
322:     }
323:   } else {
324:     SETERRQ1(PETSC_ERR_SUP, "Viewer type %s not supported for this Scotch partitioner",((PetscObject) viewer)->type_name);
325:   }
326:   return(0);
327: }

331: /*@
332:      MatPartitioningScotchSetGlobal - Set method for global partitioning.

334:   Input Parameter:
335: .  part - the partitioning context
336: .  method - MP_SCOTCH_GREED, MP_SCOTCH_GIBBS or MP_SCOTCH_GR_GI (the combination of two)
337:    Level: advanced

339: @*/
340: PetscErrorCode  MatPartitioningScotchSetGlobal(MatPartitioning part,
341:     MPScotchGlobalType global)
342: {
343:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


347:     switch (global) {
348:     case MP_SCOTCH_GREEDY:
349:     case MP_SCOTCH_GPS:
350:     case MP_SCOTCH_GR_GPS:
351:         scotch->global_method = global;
352:         break;
353:     default:
354:         SETERRQ(PETSC_ERR_SUP, "Scotch: Unknown or unsupported option");
355:     }

357:     return(0);
358: }

362: /*@
363:     MatPartitioningScotchSetCoarseLevel - Set the coarse level 
364:     
365:   Input Parameter:
366: .  part - the partitioning context
367: .  level - the coarse level in range [0.0,1.0]

369:    Level: advanced

371: @*/
372: PetscErrorCode  MatPartitioningScotchSetCoarseLevel(MatPartitioning part, PetscReal level)
373: {
374:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


378:     if (level < 0 || level > 1.0) {
379:         SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,
380:             "Scocth: level of coarsening out of range [0.0-1.0]");
381:     } else
382:         scotch->nbvtxcoarsed = (int)(part->adj->N * level);

384:     if (scotch->nbvtxcoarsed < 20)
385:         scotch->nbvtxcoarsed = 20;

387:     return(0);
388: }

392: /*@C
393:     MatPartitioningScotchSetStrategy - Set the strategy to be used by Scotch.
394:     This is an alternative way of specifying the global method, the local
395:     method, the coarse level and the multilevel option.
396:     
397:   Input Parameter:
398: .  part - the partitioning context
399: .  level - the strategy in Scotch format. Check Scotch documentation.

401:    Level: advanced

403: .seealso: MatPartitioningScotchSetGlobal(), MatPartitioningScotchSetLocal(), MatPartitioningScotchSetCoarseLevel(), MatPartitioningScotchSetMultilevel(), 
404: @*/
405: PetscErrorCode  MatPartitioningScotchSetStrategy(MatPartitioning part, char *strat)
406: {
407:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


411:     PetscStrcpy(scotch->strategy, strat);
412:     return(0);
413: }


418: /*@
419:      MatPartitioningScotchSetLocal - Set method for local partitioning.

421:   Input Parameter:
422: .  part - the partitioning context
423: .  method - MP_SCOTCH_KERNIGHAN_LIN or MP_SCOTCH_NONE

425:    Level: advanced

427: @*/
428: PetscErrorCode  MatPartitioningScotchSetLocal(MatPartitioning part, MPScotchLocalType local)
429: {
430:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


434:     switch (local) {
435:     case MP_SCOTCH_KERNIGHAN_LIN:
436:     case MP_SCOTCH_NONE:
437:         scotch->local_method = local;
438:         break;
439:     default:
440:         SETERRQ(PETSC_ERR_ARG_CORRUPT, "Scotch: Unknown or unsupported option");
441:     }

443:     return(0);
444: }

448: /*@C
449:      MatPartitioningScotchSetArch - Specify the file that describes the
450:      architecture used for mapping. The format of this file is documented in
451:      the Scotch manual.

453:   Input Parameter:
454: .  part - the partitioning context
455: .  file - the name of file
456:    Level: advanced

458:   Note:
459:   If the name is not set, then the default "archgraph.src" is used.

461: .seealso: MatPartitioningScotchSetHostList(),MatPartitioningScotchSetMapping()
462: @*/
463: PetscErrorCode  MatPartitioningScotchSetArch(MatPartitioning part, const char *filename)
464: {
465:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


469:     PetscStrcpy(scotch->arch, filename);

471:     return(0);
472: }

476: /*@C
477:      MatPartitioningScotchSetHostList - Specify host list file for mapping.

479:   Input Parameter:
480: .  part - the partitioning context
481: .  file - the name of file

483:    Level: advanced

485:   Notes:
486:   The file must consist in a list of hostnames (one per line). These hosts
487:   are the ones referred to in the architecture file (see 
488:   MatPartitioningScotchSetArch()): the first host corresponds to index 0,
489:   the second one to index 1, and so on.
490:   
491:   If the name is not set, then the default "host_list" is used.
492:   
493: .seealso: MatPartitioningScotchSetArch(), MatPartitioningScotchSetMapping()
494: @*/
495: PetscErrorCode  MatPartitioningScotchSetHostList(MatPartitioning part, const char *filename)
496: {
497:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


501:     PetscStrcpy(scotch->host_list, filename);

503:     return(0);
504: }

508: /*@
509:      MatPartitioningScotchSetMultilevel - Activates multilevel partitioning.

511:   Input Parameter:
512: .  part - the partitioning context

514:    Level: advanced

516: @*/
517: PetscErrorCode  MatPartitioningScotchSetMultilevel(MatPartitioning part)
518: {
519:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


523:     scotch->multilevel = 1;

525:     return(0);
526: }


531: /*@
532:      MatPartitioningScotchSetMapping - Activates architecture mapping for the 
533:      partitioning algorithm. Architecture mapping tries to enhance the quality
534:      of partitioning by using network topology information. 

536:   Input Parameter:
537: .  part - the partitioning context

539:    Level: advanced

541: .seealso: MatPartitioningScotchSetArch(),MatPartitioningScotchSetHostList()
542: @*/
543: PetscErrorCode  MatPartitioningScotchSetMapping(MatPartitioning part)
544: {
545:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;


549:     scotch->map = 1;

551:     return(0);
552: }

556: PetscErrorCode MatPartitioningSetFromOptions_Scotch(MatPartitioning part)
557: {
559:     PetscTruth flag;
560:     char name[PETSC_MAX_PATH_LEN];
561:     int i;
562:     PetscReal r;

564:     const char *global[] = { "greedy", "gps", "gr_gps" };
565:     const char *local[] = { "kernighan-lin", "none" };

568:     PetscOptionsHead("Set Scotch partitioning options");

570:     PetscOptionsEList("-mat_partitioning_scotch_global",
571:         "Global method to use", "MatPartitioningScotchSetGlobal", global, 3,
572:         global[0], &i, &flag);
573:     if (flag)
574:         MatPartitioningScotchSetGlobal(part, (MPScotchGlobalType)i);

576:     PetscOptionsEList("-mat_partitioning_scotch_local",
577:         "Local method to use", "MatPartitioningScotchSetLocal", local, 2,
578:         local[0], &i, &flag);
579:     if (flag)
580:         MatPartitioningScotchSetLocal(part, (MPScotchLocalType)i);

582:     PetscOptionsName("-mat_partitioning_scotch_mapping", "Use mapping",
583:         "MatPartitioningScotchSetMapping", &flag);
584:     if (flag)
585:         MatPartitioningScotchSetMapping(part);

587:     PetscOptionsString("-mat_partitioning_scotch_arch",
588:         "architecture file in scotch format", "MatPartitioningScotchSetArch",
589:         "archgraph.src", name, PETSC_MAX_PATH_LEN, &flag);
590:     if (flag)
591:         MatPartitioningScotchSetArch(part, name);

593:     PetscOptionsString("-mat_partitioning_scotch_hosts",
594:         "host list filename", "MatPartitioningScotchSetHostList",
595:         "host_list", name, PETSC_MAX_PATH_LEN, &flag);
596:     if (flag)
597:         MatPartitioningScotchSetHostList(part, name);

599:     PetscOptionsReal("-mat_partitioning_scotch_coarse_level",
600:         "coarse level", "MatPartitioningScotchSetCoarseLevel", 0, &r,
601:         &flag);
602:     if (flag)
603:         MatPartitioningScotchSetCoarseLevel(part, r);

605:     PetscOptionsName("-mat_partitioning_scotch_mul", "Use coarse level",
606:         "MatPartitioningScotchSetMultilevel", &flag);
607:     if (flag)
608:         MatPartitioningScotchSetMultilevel(part);

610:     PetscOptionsString("-mat_partitioning_scotch_strategy",
611:         "Scotch strategy string",
612:         "MatPartitioningScotchSetStrategy", "", name, PETSC_MAX_PATH_LEN,
613:         &flag);
614:     if (flag)
615:         MatPartitioningScotchSetStrategy(part, name);

617:     PetscOptionsTail();
618:     return(0);
619: }

623: PetscErrorCode MatPartitioningDestroy_Scotch(MatPartitioning part)
624: {
625:     MatPartitioning_Scotch *scotch = (MatPartitioning_Scotch *) part->data;
626:     PetscErrorCode         ierr;

629:     PetscFree(scotch->mesg_log);
630:     PetscFree(scotch);
631:     return(0);
632: }


638: /*@C
639:    MAT_PARTITIONING_SCOTCH - Creates a partitioning context via the external package SCOTCH.

641:    Collective on MPI_Comm

643:    Input Parameter:
644: .  part - the partitioning context

646:    Options Database Keys:
647: +  -mat_partitioning_scotch_global <greedy> (one of) greedy gps gr_gps
648: .  -mat_partitioning_scotch_local <kernighan-lin> (one of) kernighan-lin none
649: .  -mat_partitioning_scotch_mapping: Use mapping (MatPartitioningScotchSetMapping)
650: .  -mat_partitioning_scotch_arch <archgraph.src>: architecture file in scotch format (MatPartitioningScotchSetArch)
651: .  -mat_partitioning_scotch_hosts <host_list>: host list filename (MatPartitioningScotchSetHostList)
652: .  -mat_partitioning_scotch_coarse_level <0>: coarse level (MatPartitioningScotchSetCoarseLevel)
653: .  -mat_partitioning_scotch_mul: Use coarse level (MatPartitioningScotchSetMultilevel)
654: -  -mat_partitioning_scotch_strategy <>: Scotch strategy string (MatPartitioningScotchSetStrategy)

656:    Level: beginner

658:    Notes: See http://www.labri.fr/Perso/~pelegrin/scotch/

660: .keywords: Partitioning, create, context

662: .seealso: MatPartitioningSetType(), MatPartitioningType

664: @*/
665: PetscErrorCode  MatPartitioningCreate_Scotch(MatPartitioning part)
666: {
668:     MatPartitioning_Scotch *scotch;

671:     PetscNew(MatPartitioning_Scotch, &scotch);

673:     scotch->map = 0;
674:     scotch->global_method = MP_SCOTCH_GR_GPS;
675:     scotch->local_method = MP_SCOTCH_KERNIGHAN_LIN;
676:     PetscStrcpy(scotch->arch, "archgraph.src");
677:     scotch->nbvtxcoarsed = 200;
678:     PetscStrcpy(scotch->strategy, "");
679:     scotch->multilevel = 0;
680:     scotch->mesg_log = NULL;

682:     PetscStrcpy(scotch->host_list, "host_list");

684:     part->ops->apply = MatPartitioningApply_Scotch;
685:     part->ops->view = MatPartitioningView_Scotch;
686:     part->ops->destroy = MatPartitioningDestroy_Scotch;
687:     part->ops->setfromoptions = MatPartitioningSetFromOptions_Scotch;
688:     part->data = (void*) scotch;

690:     return(0);
691: }