Actual source code: party.c

petsc-3.3-p7 2013-05-11
  2: #include <../src/mat/impls/adj/mpi/mpiadj.h>       /*I "petscmat.h" I*/

  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: /* 
 13:    Currently using Party-1.99
 14: */
 15: EXTERN_C_BEGIN
 16: #include <party_lib.h>
 17: EXTERN_C_END

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

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

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

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

 70:   if (!flg) { /* convert regular matrix to MPIADJ */
 71:     MatConvert(matSeq,MATMPIADJ,MAT_INITIAL_MATRIX,&matAdj);
 72:   } else {
 73:     PetscObjectReference((PetscObject)matSeq);
 74:     matAdj = matSeq;
 75:   }

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

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

 90:   PetscMalloc((mat->rmap->N)*sizeof(int),&part_party);

 92:   /* redirect output to buffer */
 93: #ifdef PETSC_HAVE_UNISTD_H
 94:   fd_stdout = dup(1);
 95:   if (pipe(fd_pipe)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SYS,"Could not open pipe");
 96:   close(1);
 97:   dup2(fd_pipe[1],1);
 98:   PetscMalloc(SIZE_LOG*sizeof(char),&mesg_log);
 99: #endif

101:   /* library call */
102:   party_lib_times_start();
103:   party_lib(n,vertex_w,PETSC_NULL,PETSC_NULL,PETSC_NULL,edge_p,edge,
104:                    PETSC_NULL,p,part_party,&cutsize,redl,(char*)redm,(char*)redo,
105:                    party->global,party->local,rec,1);

107:   party_lib_times_output(1);
108:   part_info(n,vertex_w,edge_p,edge,PETSC_NULL,p,part_party,1);

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

126:   PetscMalloc((mat->rmap->N)*sizeof(PetscInt),&parttab);
127:   for (i=0;i<mat->rmap->N;i++) parttab[i] = part_party[i];

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

137:   ISCreateGeneral(((PetscObject)part)->comm,nb_locals,locals,PETSC_COPY_VALUES,partitioning);

139:   /* clean up */
140:   PetscFree(parttab);
141:   PetscFree(part_party);
142:   MatDestroy(&matSeq);
143:   MatDestroy(&matAdj);
144:   return(0);
145: }

149: PetscErrorCode MatPartitioningView_Party(MatPartitioning part,PetscViewer viewer)
150: {
151:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
152:   PetscErrorCode        ierr;
153:   PetscBool             isascii;

156:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
157:   if (isascii) {
158:     PetscViewerASCIIPrintf(viewer,"  Global method: %s\n",party->global);
159:     PetscViewerASCIIPrintf(viewer,"  Local method: %s\n",party->local);
160:     PetscViewerASCIIPrintf(viewer,"  Number of vertices for the coarse graph: %d\n",party->nbvtxcoarsed);
161:     if (party->redm) {
162:       PetscViewerASCIIPrintf(viewer,"  Using matching method for graph reduction\n");
163:     }
164:     if (party->redo) {
165:       PetscViewerASCIIPrintf(viewer,"  Using matching optimization\n");
166:     }
167:     if (party->recursive) {
168:       PetscViewerASCIIPrintf(viewer,"  Using recursive bipartitioning\n");
169:     }
170:   } else {
171:     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported for Party partitioner",((PetscObject)viewer)->type_name);
172:   }
173:   return(0);
174: }

178: /*@C
179:    MatPartitioningPartySetGlobal - Set global method for Party partitioner.

181:    Collective on MatPartitioning

183:    Input Parameters:
184: +  part - the partitioning context
185: -  method - a string representing the method

187:    Options Database:
188: .  -mat_partitioning_party_global <method> - the global method

190:    Level: advanced

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

198: .seealso: MatPartitioningPartySetLocal()
199: @*/
200: PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part,const char *global)
201: {

206:   PetscTryMethod(part,"MatPartitioningPartySetGlobal_C",(MatPartitioning,const char*),(part,global));
207:   return(0);
208: }

210: EXTERN_C_BEGIN
213: PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part,const char* global)
214: {
215:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
216:   PetscErrorCode        ierr;

219:   PetscStrncpy(party->global,global,15);
220:   return(0);
221: }
222: EXTERN_C_END

226: /*@C
227:    MatPartitioningPartySetLocal - Set local method for Party partitioner.

229:    Collective on MatPartitioning

231:    Input Parameters:
232: +  part - the partitioning context
233: -  method - a string representing the method

235:    Options Database:
236: .  -mat_partitioning_party_local <method> - the local method

238:    Level: advanced

240:    Notes:
241:    The method may be one of MP_PARTY_HELPFUL_SETS, MP_PARTY_KERNIGHAN_LIN, or
242:    MP_PARTY_NONE. Check the Party Library Users Manual for details.

244: .seealso: MatPartitioningPartySetGlobal()
245: @*/
246: PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part,const char *local)
247: {

252:   PetscTryMethod(part,"MatPartitioningPartySetLocal_C",(MatPartitioning,const char*),(part,local));
253:   return(0);
254: }

256: EXTERN_C_BEGIN
259: PetscErrorCode MatPartitioningPartySetLocal_Party(MatPartitioning part,const char* local)

261: {
262:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
263:   PetscErrorCode        ierr;

266:   PetscStrncpy(party->local,local,15);
267:   return(0);
268: }
269: EXTERN_C_END

273: /*@
274:    MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
275:    Party partitioner.
276:  
277:    Collective on MatPartitioning

279:    Input Parameters:
280: +  part - the partitioning context
281: -  level - the coarse level in range [0.0,1.0]

283:    Options Database:
284: .  -mat_partitioning_party_coarse <l> - Coarse level

286:    Level: advanced
287: @*/
288: PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part,PetscReal level)
289: {

295:   PetscTryMethod(part,"MatPartitioningPartySetCoarseLevel_C",(MatPartitioning,PetscReal),(part,level));
296:   return(0);
297: }

299: EXTERN_C_BEGIN
302: PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part,PetscReal level)
303: {
304:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

307:   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]");
308:   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
309:   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
310:   return(0);
311: }
312: EXTERN_C_END

316: /*@
317:    MatPartitioningPartySetMatchOptimization - Activate matching optimization for
318:    graph reduction.
319:  
320:    Collective on MatPartitioning

322:    Input Parameters:
323: +  part - the partitioning context
324: -  opt - boolean flag

326:    Options Database:
327: .  -mat_partitioning_party_match_optimization - Matching optimization on/off

329:    Level: advanced
330: @*/
331: PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part,PetscBool opt)
332: {

338:   PetscTryMethod(part,"MatPartitioningPartySetMatchOptimization_C",(MatPartitioning,PetscBool),(part,opt));
339:   return(0);
340: }

342: EXTERN_C_BEGIN
345: PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part,PetscBool opt)
346: {
347:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

350:   party->redo = opt;
351:   return(0);
352: }
353: EXTERN_C_END

357: /*@
358:    MatPartitioningPartySetBipart - Activate or deactivate recursive bisection.
359:  
360:    Collective on MatPartitioning

362:    Input Parameters:
363: +  part - the partitioning context
364: -  bp - boolean flag

366:    Options Database:
367: -  -mat_partitioning_party_bipart - Bipartitioning option on/off

369:    Level: advanced
370: @*/
371: PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part,PetscBool bp)
372: {

378:   PetscTryMethod(part,"MatPartitioningPartySetBipart_C",(MatPartitioning,PetscBool),(part,bp));
379:   return(0);
380: }

382: EXTERN_C_BEGIN
385: PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part,PetscBool bp)
386: {
387:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

390:   party->recursive = bp;
391:   return(0);
392: }
393: EXTERN_C_END

397: PetscErrorCode MatPartitioningSetFromOptions_Party(MatPartitioning part)
398: {
399:   PetscErrorCode        ierr;
400:   PetscBool             flag;
401:   char                  value[256];
402:   PetscReal             r;
403:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

406:   PetscOptionsHead("Set Party partitioning options");
407:     PetscOptionsString("-mat_partitioning_party_global","Global method","MatPartitioningPartySetGlobal",party->global,value,256,&flag);
408:     if (flag) { MatPartitioningPartySetGlobal(part,value); }
409:     PetscOptionsString("-mat_partitioning_party_local","Local method","MatPartitioningPartySetLocal",party->local,value,256,&flag);
410:     if (flag) { MatPartitioningPartySetLocal(part,value); }
411:     PetscOptionsReal("-mat_partitioning_party_coarse","Coarse level","MatPartitioningPartySetCoarseLevel",0.0,&r,&flag);
412:     if (flag) { MatPartitioningPartySetCoarseLevel(part,r); }
413:     PetscOptionsBool("-mat_partitioning_party_match_optimization","Matching optimization on/off","MatPartitioningPartySetMatchOptimization",party->redo,&party->redo,PETSC_NULL);
414:     PetscOptionsBool("-mat_partitioning_party_bipart","Bipartitioning on/off","MatPartitioningPartySetBipart",party->recursive,&party->recursive,PETSC_NULL);
415:     PetscOptionsBool("-mat_partitioning_party_verbose","Show library output","",party->verbose,&party->verbose,PETSC_NULL);
416:   PetscOptionsTail();
417:   return(0);
418: }

422: PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
423: {
424:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
425:   PetscErrorCode        ierr;

428:   PetscFree(party);
429:   /* clear composed functions */
430:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetGlobal_C","",PETSC_NULL);
431:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetLocal_C","",PETSC_NULL);
432:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C","",PETSC_NULL);
433:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C","",PETSC_NULL);
434:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetBipart_C","",PETSC_NULL);
435:   return(0);
436: }

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

441:    Level: beginner

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

445: .keywords: Partitioning, create, context

447: .seealso: MatPartitioningSetType(), MatPartitioningType

449: M*/

451: EXTERN_C_BEGIN
454: PetscErrorCode  MatPartitioningCreate_Party(MatPartitioning part)
455: {
456:   PetscErrorCode        ierr;
457:   MatPartitioning_Party *party;

460:   PetscNewLog(part,MatPartitioning_Party,&party);
461:   part->data = (void*)party;

463:   PetscStrcpy(party->global,"gcf,gbf");
464:   PetscStrcpy(party->local,"kl");
465:   party->redm         = PETSC_TRUE;
466:   party->redo         = PETSC_TRUE;
467:   party->recursive    = PETSC_TRUE;
468:   party->verbose      = PETSC_FALSE;
469:   party->nbvtxcoarsed = 200;

471:   part->ops->apply          = MatPartitioningApply_Party;
472:   part->ops->view           = MatPartitioningView_Party;
473:   part->ops->destroy        = MatPartitioningDestroy_Party;
474:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;

476:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetGlobal_C","MatPartitioningPartySetGlobal_Party",MatPartitioningPartySetGlobal_Party);
477:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetLocal_C","MatPartitioningPartySetLocal_Party",MatPartitioningPartySetLocal_Party);
478:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C","MatPartitioningPartySetCoarseLevel_Party",MatPartitioningPartySetCoarseLevel_Party);
479:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C","MatPartitioningPartySetMatchOptimization_Party",MatPartitioningPartySetMatchOptimization_Party);
480:   PetscObjectComposeFunctionDynamic((PetscObject)part,"MatPartitioningPartySetBipart_C","MatPartitioningPartySetBipart_Party",MatPartitioningPartySetBipart_Party);
481:   return(0);
482: }
483: EXTERN_C_END