Actual source code: party.c


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

  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 */

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

 45:   MPI_Comm_size(PetscObjectComm((PetscObject)mat),&size);
 46:   MPI_Comm_rank(PetscObjectComm((PetscObject)mat),&rank);
 47:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
 48:   if (size>1) {
 49:     if (flg) {
 50:       MatMPIAdjToSeq(mat,&matSeq);
 51:      } else {
 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:       MatCreateSubMatrices(mat,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&A);
 57:       ISDestroy(&isrow);
 58:       ISDestroy(&iscol);
 59:       matSeq = *A;
 60:       PetscFree(A);
 61:      }
 62:   } else {
 63:     PetscObjectReference((PetscObject)mat);
 64:     matSeq = mat;
 65:   }

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

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

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

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

 89:   /* redirect output to buffer */
 90: #if defined(PETSC_HAVE_UNISTD_H)
 91:   fd_stdout = dup(1);
 93:   close(1);
 94:   dup2(fd_pipe[1],1);
 95:   PetscMalloc1(SIZE_LOG,&mesg_log);
 96: #endif

 98:   /* library call */
 99:   party_lib_times_start();
100:   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);

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

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

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

126:   /* creation of the index set */
127:   nb_locals = mat->rmap->n;
128:   locals    = parttab + mat->rmap->rstart;

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

132:   /* clean up */
133:   PetscFree(parttab);
134:   PetscFree(part_party);
135:   MatDestroy(&matSeq);
136:   MatDestroy(&matAdj);
137:   return 0;
138: }

140: PetscErrorCode MatPartitioningView_Party(MatPartitioning part,PetscViewer viewer)
141: {
142:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;
143:   PetscBool             isascii;

145:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
146:   if (isascii) {
147:     PetscViewerASCIIPrintf(viewer,"  Global method: %s\n",party->global);
148:     PetscViewerASCIIPrintf(viewer,"  Local method: %s\n",party->local);
149:     PetscViewerASCIIPrintf(viewer,"  Number of vertices for the coarse graph: %d\n",party->nbvtxcoarsed);
150:     if (party->redm) PetscViewerASCIIPrintf(viewer,"  Using matching method for graph reduction\n");
151:     if (party->redo) PetscViewerASCIIPrintf(viewer,"  Using matching optimization\n");
152:     if (party->recursive) PetscViewerASCIIPrintf(viewer,"  Using recursive bipartitioning\n");
153:   }
154:   return 0;
155: }

157: /*@C
158:    MatPartitioningPartySetGlobal - Set global method for Party partitioner.

160:    Collective on MatPartitioning

162:    Input Parameters:
163: +  part - the partitioning context
164: -  method - a string representing the method

166:    Options Database:
167: .  -mat_partitioning_party_global <method> - the global method

169:    Level: advanced

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

177: .seealso: MatPartitioningPartySetLocal()
178: @*/
179: PetscErrorCode MatPartitioningPartySetGlobal(MatPartitioning part,const char *global)
180: {
182:   PetscTryMethod(part,"MatPartitioningPartySetGlobal_C",(MatPartitioning,const char*),(part,global));
183:   return 0;
184: }

186: PetscErrorCode MatPartitioningPartySetGlobal_Party(MatPartitioning part,const char *global)
187: {
188:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

190:   PetscStrncpy(party->global,global,15);
191:   return 0;
192: }

194: /*@C
195:    MatPartitioningPartySetLocal - Set local method for Party partitioner.

197:    Collective on MatPartitioning

199:    Input Parameters:
200: +  part - the partitioning context
201: -  method - a string representing the method

203:    Options Database:
204: .  -mat_partitioning_party_local <method> - the local method

206:    Level: advanced

208:    Notes:
209:    The method may be one of MP_PARTY_HELPFUL_SETS, MP_PARTY_KERNIGHAN_LIN, or
210:    MP_PARTY_NONE. Check the Party Library Users Manual for details.

212: .seealso: MatPartitioningPartySetGlobal()
213: @*/
214: PetscErrorCode MatPartitioningPartySetLocal(MatPartitioning part,const char *local)
215: {
217:   PetscTryMethod(part,"MatPartitioningPartySetLocal_C",(MatPartitioning,const char*),(part,local));
218:   return 0;
219: }

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

223: {
224:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

226:   PetscStrncpy(party->local,local,15);
227:   return 0;
228: }

230: /*@
231:    MatPartitioningPartySetCoarseLevel - Set the coarse level parameter for the
232:    Party partitioner.

234:    Collective on MatPartitioning

236:    Input Parameters:
237: +  part - the partitioning context
238: -  level - the coarse level in range [0.0,1.0]

240:    Options Database:
241: .  -mat_partitioning_party_coarse <l> - Coarse level

243:    Level: advanced
244: @*/
245: PetscErrorCode MatPartitioningPartySetCoarseLevel(MatPartitioning part,PetscReal level)
246: {
249:   PetscTryMethod(part,"MatPartitioningPartySetCoarseLevel_C",(MatPartitioning,PetscReal),(part,level));
250:   return 0;
251: }

253: PetscErrorCode MatPartitioningPartySetCoarseLevel_Party(MatPartitioning part,PetscReal level)
254: {
255:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

258:   party->nbvtxcoarsed = (PetscInt)(part->adj->cmap->N * level);
259:   if (party->nbvtxcoarsed < 20) party->nbvtxcoarsed = 20;
260:   return 0;
261: }

263: /*@
264:    MatPartitioningPartySetMatchOptimization - Activate matching optimization for
265:    graph reduction.

267:    Collective on MatPartitioning

269:    Input Parameters:
270: +  part - the partitioning context
271: -  opt - boolean flag

273:    Options Database:
274: .  -mat_partitioning_party_match_optimization - Matching optimization on/off

276:    Level: advanced
277: @*/
278: PetscErrorCode MatPartitioningPartySetMatchOptimization(MatPartitioning part,PetscBool opt)
279: {
282:   PetscTryMethod(part,"MatPartitioningPartySetMatchOptimization_C",(MatPartitioning,PetscBool),(part,opt));
283:   return 0;
284: }

286: PetscErrorCode MatPartitioningPartySetMatchOptimization_Party(MatPartitioning part,PetscBool opt)
287: {
288:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

290:   party->redo = opt;
291:   return 0;
292: }

294: /*@
295:    MatPartitioningPartySetBipart - Activate or deactivate recursive bisection.

297:    Collective on MatPartitioning

299:    Input Parameters:
300: +  part - the partitioning context
301: -  bp - boolean flag

303:    Options Database:
304: -  -mat_partitioning_party_bipart - Bipartitioning option on/off

306:    Level: advanced
307: @*/
308: PetscErrorCode MatPartitioningPartySetBipart(MatPartitioning part,PetscBool bp)
309: {
312:   PetscTryMethod(part,"MatPartitioningPartySetBipart_C",(MatPartitioning,PetscBool),(part,bp));
313:   return 0;
314: }

316: PetscErrorCode MatPartitioningPartySetBipart_Party(MatPartitioning part,PetscBool bp)
317: {
318:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

320:   party->recursive = bp;
321:   return 0;
322: }

324: PetscErrorCode MatPartitioningSetFromOptions_Party(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
325: {
326:   PetscBool             flag;
327:   char                  value[256];
328:   PetscReal             r;
329:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

331:   PetscOptionsHead(PetscOptionsObject,"Set Party partitioning options");
332:   PetscOptionsString("-mat_partitioning_party_global","Global method","MatPartitioningPartySetGlobal",party->global,value,sizeof(value),&flag);
333:   if (flag) MatPartitioningPartySetGlobal(part,value);
334:   PetscOptionsString("-mat_partitioning_party_local","Local method","MatPartitioningPartySetLocal",party->local,value,sizeof(value),&flag);
335:   if (flag) MatPartitioningPartySetLocal(part,value);
336:   PetscOptionsReal("-mat_partitioning_party_coarse","Coarse level","MatPartitioningPartySetCoarseLevel",0.0,&r,&flag);
337:   if (flag) MatPartitioningPartySetCoarseLevel(part,r);
338:   PetscOptionsBool("-mat_partitioning_party_match_optimization","Matching optimization on/off","MatPartitioningPartySetMatchOptimization",party->redo,&party->redo,NULL);
339:   PetscOptionsBool("-mat_partitioning_party_bipart","Bipartitioning on/off","MatPartitioningPartySetBipart",party->recursive,&party->recursive,NULL);
340:   PetscOptionsBool("-mat_partitioning_party_verbose","Show library output","",party->verbose,&party->verbose,NULL);
341:   PetscOptionsTail();
342:   return 0;
343: }

345: PetscErrorCode MatPartitioningDestroy_Party(MatPartitioning part)
346: {
347:   MatPartitioning_Party *party = (MatPartitioning_Party*)part->data;

349:   PetscFree(party);
350:   /* clear composed functions */
351:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetGlobal_C",NULL);
352:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetLocal_C",NULL);
353:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C",NULL);
354:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C",NULL);
355:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetBipart_C",NULL);
356:   return 0;
357: }

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

362:    Level: beginner

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

367:     Does not support using MatPartitioningSetUseEdgeWeights()

369: .seealso: MatPartitioningSetType(), MatPartitioningType

371: M*/

373: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Party(MatPartitioning part)
374: {
375:   MatPartitioning_Party *party;

377:   PetscNewLog(part,&party);
378:   part->data = (void*)party;

380:   PetscStrcpy(party->global,"gcf,gbf");
381:   PetscStrcpy(party->local,"kl");

383:   party->redm         = PETSC_TRUE;
384:   party->redo         = PETSC_TRUE;
385:   party->recursive    = PETSC_TRUE;
386:   party->verbose      = PETSC_FALSE;
387:   party->nbvtxcoarsed = 200;

389:   part->ops->apply          = MatPartitioningApply_Party;
390:   part->ops->view           = MatPartitioningView_Party;
391:   part->ops->destroy        = MatPartitioningDestroy_Party;
392:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Party;

394:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetGlobal_C",MatPartitioningPartySetGlobal_Party);
395:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetLocal_C",MatPartitioningPartySetLocal_Party);
396:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetCoarseLevel_C",MatPartitioningPartySetCoarseLevel_Party);
397:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetMatchOptimization_C",MatPartitioningPartySetMatchOptimization_Party);
398:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPartySetBipart_C",MatPartitioningPartySetBipart_Party);
399:   return 0;
400: }