Actual source code: party.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: #include <../src/mat/impls/adj/mpi/mpiadj.h>       /*I "petscmat.h" I*/

  4: #if defined(PETSC_HAVE_UNISTD_H)
  5: #include <unistd.h>
  6: #endif

  8: /*
  9:    Currently using Party-1.99
 10: */
 11: EXTERN_C_BEGIN
 12: #include <party_lib.h>
 13: EXTERN_C_END

 15: typedef struct {
 16:   PetscBool redm;
 17:   PetscBool redo;
 18:   PetscBool recursive;
 19:   PetscBool verbose;
 20:   char      global[15];         /* global method */
 21:   char      local[15];          /* local method */
 22:   PetscInt  nbvtxcoarsed;       /* number of vertices for the coarse graph */
 23: } MatPartitioning_Party;

 25: #define SIZE_LOG 10000          /* size of buffer for mesg_log */

 29: static PetscErrorCode MatPartitioningApply_Party(MatPartitioning part,IS *partitioning)
 30: {
 31:   PetscErrorCode        ierr;
 32:   PetscInt              i,*parttab,*locals,nb_locals,M,N;
 33:   PetscMPIInt           size,rank;
 34:   Mat                   mat = part->adj,matAdj,matSeq,*A;
 35:   Mat_MPIAdj            *adj;
 36:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
 37:   PetscBool             flg;
 38:   IS                    isrow, iscol;
 39:   int                   n,*edge_p,*edge,*vertex_w,p,*part_party,cutsize,redl,rec;
 40:   const char            *redm,*redo;
 41:   char                  *mesg_log;
 42: #if defined(PETSC_HAVE_UNISTD_H)
 43:   int                   fd_stdout,fd_pipe[2],count,err;
 44: #endif

 47:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
 48:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
 49:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
 50:   if (size>1) {
 51:     if (flg) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Distributed matrix format MPIAdj is not supported for sequential partitioners");
 52:     PetscInfo(part,"Converting distributed matrix to sequential: this could be a performance loss\n");
 53:     MatGetSize(mat,&M,&N);
 54:     ISCreateStride(PETSC_COMM_SELF,M,0,1,&isrow);
 55:     ISCreateStride(PETSC_COMM_SELF,N,0,1,&iscol);
 56:     MatGetSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&A);
 57:     ISDestroy(&isrow);
 58:     ISDestroy(&iscol);
 59:     matSeq = *A;
 60:     PetscFree(A);
 61:   } else {
 62:     PetscObjectReference((PetscObject)mat);
 63:     matSeq = mat;
 64:   }

 66:   if (!flg) { /* convert regular matrix to MPIADJ */
 67:     MatConvert(matSeq,MATMPIADJ,MAT_INITIAL_MATRIX,&matAdj);
 68:   } else {
 69:     PetscObjectReference((PetscObject)matSeq);
 70:     matAdj = matSeq;
 71:   }

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

 75:   /* arguments for Party library */
 76:   n        = mat->rmap->N;             /* number of vertices in full graph */
 77:   edge_p   = adj->i;                   /* start of edge list for each vertex */
 78:   edge     = adj->j;                   /* edge list data */
 79:   vertex_w = part->vertex_weights;     /* weights for all vertices */
 80:   p        = part->n;                  /* number of parts to create */
 81:   redl     = party->nbvtxcoarsed;      /* how many vertices to coarsen down to? */
 82:   rec      = party->recursive ? 1 : 0; /* recursive bisection */
 83:   redm     = party->redm ? "lam" : ""; /* matching method */
 84:   redo     = party->redo ? "w3" : "";  /* matching optimization method */

 86:   PetscMalloc1(mat->rmap->N,&part_party);

 88:   /* redirect output to buffer */
 89: #if defined(PETSC_HAVE_UNISTD_H)
 90:   fd_stdout = dup(1);
 91:   if (pipe(fd_pipe)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Could not open pipe");
 92:   close(1);
 93:   dup2(fd_pipe[1],1);
 94:   PetscMalloc1(SIZE_LOG,&mesg_log);
 95: #endif

 97:   /* library call */
 98:   party_lib_times_start();
 99:   party_lib(n,vertex_w,NULL,NULL,NULL,edge_p,edge,NULL,p,part_party,&cutsize,redl,(char*)redm,(char*)redo,party->global,party->local,rec,1);

101:   party_lib_times_output(1);
102:   part_info(n,vertex_w,edge_p,edge,NULL,p,part_party,1);

104: #if defined(PETSC_HAVE_UNISTD_H)
105:   err = fflush(stdout);
106:   if (err) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"fflush() failed on stdout");
107:   count = read(fd_pipe[0],mesg_log,(SIZE_LOG-1)*sizeof(char));
108:   if (count<0) count = 0;
109:   mesg_log[count] = 0;
110:   close(1);
111:   dup2(fd_stdout,1);
112:   close(fd_stdout);
113:   close(fd_pipe[0]);
114:   close(fd_pipe[1]);
115:   if (party->verbose) {
116:     PetscPrintf(PetscObjectComm((PetscObject)mat),mesg_log);
117:   }
118:   PetscFree(mesg_log);
119: #endif
120:   if (ierr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Party failed");

122:   PetscMalloc1(mat->rmap->N,&parttab);
123:   for (i=0; i<mat->rmap->N; i++) parttab[i] = part_party[i];

125:   /* creation of the index set */
126:   nb_locals = mat->rmap->N / size;
127:   locals    = parttab + rank*nb_locals;
128:   if (rank < mat->rmap->N % size) {
129:     nb_locals++;
130:     locals += rank;
131:   } else locals += mat->rmap->N % size;

133:   ISCreateGeneral(PetscObjectComm((PetscObject)part),nb_locals,locals,PETSC_COPY_VALUES,partitioning);

135:   /* clean up */
136:   PetscFree(parttab);
137:   PetscFree(part_party);
138:   MatDestroy(&matSeq);
139:   MatDestroy(&matAdj);
140:   return(0);
141: }

145: PetscErrorCode MatPartitioningView_Party(MatPartitioning part,PetscViewer viewer)
146: {
147:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
148:   PetscErrorCode        ierr;
149:   PetscBool             isascii;

152:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
153:   if (isascii) {
154:     PetscViewerASCIIPrintf(viewer,"  Global method: %s\n",party->global);
155:     PetscViewerASCIIPrintf(viewer,"  Local method: %s\n",party->local);
156:     PetscViewerASCIIPrintf(viewer,"  Number of vertices for the coarse graph: %d\n",party->nbvtxcoarsed);
157:     if (party->redm) {
158:       PetscViewerASCIIPrintf(viewer,"  Using matching method for graph reduction\n");
159:     }
160:     if (party->redo) {
161:       PetscViewerASCIIPrintf(viewer,"  Using matching optimization\n");
162:     }
163:     if (party->recursive) {
164:       PetscViewerASCIIPrintf(viewer,"  Using recursive bipartitioning\n");
165:     }
166:   }
167:   return(0);
168: }

172: /*@C
173:    MatPartitioningPartySetGlobal - Set global method for Party partitioner.

175:    Collective on MatPartitioning

177:    Input Parameters:
178: +  part - the partitioning context
179: -  method - a string representing the method

181:    Options Database:
182: .  -mat_partitioning_party_global <method> - the global method

184:    Level: advanced

186:    Notes:
187:    The method may be one of MP_PARTY_OPT, MP_PARTY_LIN, MP_PARTY_SCA,
188:    MP_PARTY_RAN, MP_PARTY_GBF, MP_PARTY_GCF, MP_PARTY_BUB or MP_PARTY_DEF, or
189:    alternatively a string describing the method. Two or more methods can be
190:    combined like "gbf,gcf". Check the Party Library Users Manual for details.

192: .seealso: MatPartitioningPartySetLocal()
193: @*/
194: PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part,const char *global)
195: {

200:   PetscTryMethod(part,"MatPartitioningPartySetGlobal_C",(MatPartitioning,const char*),(part,global));
201:   return(0);
202: }

206: PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part,const char *global)
207: {
208:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
209:   PetscErrorCode        ierr;

212:   PetscStrncpy(party->global,global,15);
213:   return(0);
214: }

218: /*@C
219:    MatPartitioningPartySetLocal - Set local method for Party partitioner.

221:    Collective on MatPartitioning

223:    Input Parameters:
224: +  part - the partitioning context
225: -  method - a string representing the method

227:    Options Database:
228: .  -mat_partitioning_party_local <method> - the local method

230:    Level: advanced

232:    Notes:
233:    The method may be one of MP_PARTY_HELPFUL_SETS, MP_PARTY_KERNIGHAN_LIN, or
234:    MP_PARTY_NONE. Check the Party Library Users Manual for details.

236: .seealso: MatPartitioningPartySetGlobal()
237: @*/
238: PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part,const char *local)
239: {

244:   PetscTryMethod(part,"MatPartitioningPartySetLocal_C",(MatPartitioning,const char*),(part,local));
245:   return(0);
246: }

250: PetscErrorCode MatPartitioningPartySetLocal_Party(MatPartitioning part,const char *local)

252: {
253:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
254:   PetscErrorCode        ierr;

257:   PetscStrncpy(party->local,local,15);
258:   return(0);
259: }

263: /*@
264:    MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
265:    Party partitioner.

267:    Collective on MatPartitioning

269:    Input Parameters:
270: +  part - the partitioning context
271: -  level - the coarse level in range [0.0,1.0]

273:    Options Database:
274: .  -mat_partitioning_party_coarse <l> - Coarse level

276:    Level: advanced
277: @*/
278: PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part,PetscReal level)
279: {

285:   PetscTryMethod(part,"MatPartitioningPartySetCoarseLevel_C",(MatPartitioning,PetscReal),(part,level));
286:   return(0);
287: }

291: PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part,PetscReal level)
292: {
293:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

296:   if (level<0.0 || level>1.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Party: level of coarsening out of range [0.0-1.0]");
297:   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
298:   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
299:   return(0);
300: }

304: /*@
305:    MatPartitioningPartySetMatchOptimization - Activate matching optimization for
306:    graph reduction.

308:    Collective on MatPartitioning

310:    Input Parameters:
311: +  part - the partitioning context
312: -  opt - boolean flag

314:    Options Database:
315: .  -mat_partitioning_party_match_optimization - Matching optimization on/off

317:    Level: advanced
318: @*/
319: PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part,PetscBool opt)
320: {

326:   PetscTryMethod(part,"MatPartitioningPartySetMatchOptimization_C",(MatPartitioning,PetscBool),(part,opt));
327:   return(0);
328: }

332: PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part,PetscBool opt)
333: {
334:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

337:   party->redo = opt;
338:   return(0);
339: }

343: /*@
344:    MatPartitioningPartySetBipart - Activate or deactivate recursive bisection.

346:    Collective on MatPartitioning

348:    Input Parameters:
349: +  part - the partitioning context
350: -  bp - boolean flag

352:    Options Database:
353: -  -mat_partitioning_party_bipart - Bipartitioning option on/off

355:    Level: advanced
356: @*/
357: PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part,PetscBool bp)
358: {

364:   PetscTryMethod(part,"MatPartitioningPartySetBipart_C",(MatPartitioning,PetscBool),(part,bp));
365:   return(0);
366: }

370: PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part,PetscBool bp)
371: {
372:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

375:   party->recursive = bp;
376:   return(0);
377: }

381: PetscErrorCode MatPartitioningSetFromOptions_Party(PetscOptions *PetscOptionsObject,MatPartitioning part)
382: {
383:   PetscErrorCode        ierr;
384:   PetscBool             flag;
385:   char                  value[256];
386:   PetscReal             r;
387:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

390:   PetscOptionsHead(PetscOptionsObject,"Set Party partitioning options");
391:   PetscOptionsString("-mat_partitioning_party_global","Global method","MatPartitioningPartySetGlobal",party->global,value,256,&flag);
392:   if (flag) { MatPartitioningPartySetGlobal(part,value); }
393:   PetscOptionsString("-mat_partitioning_party_local","Local method","MatPartitioningPartySetLocal",party->local,value,256,&flag);
394:   if (flag) { MatPartitioningPartySetLocal(part,value); }
395:   PetscOptionsReal("-mat_partitioning_party_coarse","Coarse level","MatPartitioningPartySetCoarseLevel",0.0,&r,&flag);
396:   if (flag) { MatPartitioningPartySetCoarseLevel(part,r); }
397:   PetscOptionsBool("-mat_partitioning_party_match_optimization","Matching optimization on/off","MatPartitioningPartySetMatchOptimization",party->redo,&party->redo,NULL);
398:   PetscOptionsBool("-mat_partitioning_party_bipart","Bipartitioning on/off","MatPartitioningPartySetBipart",party->recursive,&party->recursive,NULL);
399:   PetscOptionsBool("-mat_partitioning_party_verbose","Show library output","",party->verbose,&party->verbose,NULL);
400:   PetscOptionsTail();
401:   return(0);
402: }

406: PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
407: {
408:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
409:   PetscErrorCode        ierr;

412:   PetscFree(party);
413:   /* clear composed functions */
414:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetGlobal_C",NULL);
415:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetLocal_C",NULL);
416:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C",NULL);
417:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C",NULL);
418:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetBipart_C",NULL);
419:   return(0);
420: }

422: /*MC
423:    MATPARTITIONINGPARTY - Creates a partitioning context via the external package Party.

425:    Level: beginner

427:    Notes: See http://wwwcs.upb.de/fachbereich/AG/monien/RESEARCH/PART/party.html

429: .keywords: Partitioning, create, context

431: .seealso: MatPartitioningSetType(), MatPartitioningType

433: M*/

437: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning part)
438: {
439:   PetscErrorCode        ierr;
440:   MatPartitioning_Party *party;

443:   PetscNewLog(part,&party);
444:   part->data = (void*)party;

446:   PetscStrcpy(party->global,"gcf,gbf");
447:   PetscStrcpy(party->local,"kl");

449:   party->redm         = PETSC_TRUE;
450:   party->redo         = PETSC_TRUE;
451:   party->recursive    = PETSC_TRUE;
452:   party->verbose      = PETSC_FALSE;
453:   party->nbvtxcoarsed = 200;

455:   part->ops->apply          = MatPartitioningApply_Party;
456:   part->ops->view           = MatPartitioningView_Party;
457:   part->ops->destroy        = MatPartitioningDestroy_Party;
458:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;

460:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetGlobal_C",MatPartitioningPartySetGlobal_Party);
461:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetLocal_C",MatPartitioningPartySetLocal_Party);
462:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C",MatPartitioningPartySetCoarseLevel_Party);
463:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C",MatPartitioningPartySetMatchOptimization_Party);
464:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetBipart_C",MatPartitioningPartySetBipart_Party);
465:   return(0);
466: }