Actual source code: scotch.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  2: #include <../src/mat/impls/adj/mpi/mpiadj.h>       /*I "petscmat.h" I*/

  4: EXTERN_C_BEGIN
  5: #include <ptscotch.h>
  6: EXTERN_C_END

  8: typedef struct {
  9:   double     imbalance;
 10:   SCOTCH_Num strategy;
 11: } MatPartitioning_PTScotch;

 15: /*@
 16:    MatPartitioningPTScotchSetImbalance - Sets the value of the load imbalance
 17:    ratio to be used during strategy selection.

 19:    Collective on MatPartitioning

 21:    Input Parameters:
 22: +  part - the partitioning context
 23: -  imb  - the load imbalance ratio

 25:    Options Database:
 26: .  -mat_partitioning_ptscotch_imbalance <imb>

 28:    Note:
 29:    Must be in the range [0,1]. The default value is 0.01.

 31:    Level: advanced

 33: .seealso: MatPartitioningPTScotchSetStrategy(), MatPartitioningPTScotchGetImbalance()
 34: @*/
 35: PetscErrorCode MatPartitioningPTScotchSetImbalance(MatPartitioning part,PetscReal imb)
 36: {

 42:   PetscTryMethod(part,"MatPartitioningPTScotchSetImbalance_C",(MatPartitioning,PetscReal),(part,imb));
 43:   return(0);
 44: }

 48: PetscErrorCode MatPartitioningPTScotchSetImbalance_PTScotch(MatPartitioning part,PetscReal imb)
 49: {
 50:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

 53:   if (imb==PETSC_DEFAULT) scotch->imbalance = 0.01;
 54:   else {
 55:     if (imb<0.0 || imb>1.0) SETERRQ(PetscObjectComm((PetscObject)part),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of imb. Must be in range [0,1]");
 56:     scotch->imbalance = (double)imb;
 57:   }
 58:   return(0);
 59: }

 63: /*@
 64:    MatPartitioningPTScotchGetImbalance - Gets the value of the load imbalance
 65:    ratio used during strategy selection.

 67:    Not Collective

 69:    Input Parameter:
 70: .  part - the partitioning context

 72:    Output Parameter:
 73: .  imb  - the load imbalance ratio

 75:    Level: advanced

 77: .seealso: MatPartitioningPTScotchSetImbalance()
 78: @*/
 79: PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning part,PetscReal *imb)
 80: {

 86:   PetscUseMethod(part,"MatPartitioningPTScotchGetImbalance_C",(MatPartitioning,PetscReal*),(part,imb));
 87:   return(0);
 88: }

 92: PetscErrorCode MatPartitioningPTScotchGetImbalance_PTScotch(MatPartitioning part,PetscReal *imb)
 93: {
 94:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

 97:   *imb = scotch->imbalance;
 98:   return(0);
 99: }

103: /*@
104:    MatPartitioningPTScotchSetStrategy - Sets the strategy to be used in PTScotch.

106:    Collective on MatPartitioning

108:    Input Parameters:
109: +  part - the partitioning context
110: -  strategy - the strategy, one of
111: .vb
112:      MP_PTSCOTCH_QUALITY     - Prioritize quality over speed
113:      MP_PTSCOTCH_SPEED       - Prioritize speed over quality
114:      MP_PTSCOTCH_BALANCE     - Enforce load balance
115:      MP_PTSCOTCH_SAFETY      - Avoid methods that may fail
116:      MP_PTSCOTCH_SCALABILITY - Favor scalability as much as possible
117: .ve

119:    Options Database:
120: .  -mat_partitioning_ptscotch_strategy [quality,speed,balance,safety,scalability] - strategy

122:    Level: advanced

124:    Notes:
125:    The default is MP_SCOTCH_QUALITY. See the PTScotch documentation for more information.

127: .seealso: MatPartitioningPTScotchSetImbalance(), MatPartitioningPTScotchGetStrategy()
128: @*/
129: PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning part,MPPTScotchStrategyType strategy)
130: {

136:   PetscTryMethod(part,"MatPartitioningPTScotchSetStrategy_C",(MatPartitioning,MPPTScotchStrategyType),(part,strategy));
137:   return(0);
138: }

142: PetscErrorCode MatPartitioningPTScotchSetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType strategy)
143: {
144:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

147:   switch (strategy) {
148:   case MP_PTSCOTCH_QUALITY:     scotch->strategy = SCOTCH_STRATQUALITY; break;
149:   case MP_PTSCOTCH_SPEED:       scotch->strategy = SCOTCH_STRATSPEED; break;
150:   case MP_PTSCOTCH_BALANCE:     scotch->strategy = SCOTCH_STRATBALANCE; break;
151:   case MP_PTSCOTCH_SAFETY:      scotch->strategy = SCOTCH_STRATSAFETY; break;
152:   case MP_PTSCOTCH_SCALABILITY: scotch->strategy = SCOTCH_STRATSCALABILITY; break;
153:   }
154:   return(0);
155: }

159: /*@
160:    MatPartitioningPTScotchGetStrategy - Gets the strategy used in PTScotch.

162:    Not Collective

164:    Input Parameter:
165: .  part - the partitioning context

167:    Output Parameter:
168: .  strategy - the strategy

170:    Level: advanced

172: .seealso: MatPartitioningPTScotchSetStrategy()
173: @*/
174: PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning part,MPPTScotchStrategyType *strategy)
175: {

181:   PetscUseMethod(part,"MatPartitioningPTScotchGetStrategy_C",(MatPartitioning,MPPTScotchStrategyType*),(part,strategy));
182:   return(0);
183: }

187: PetscErrorCode MatPartitioningPTScotchGetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType *strategy)
188: {
189:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

192:   switch (scotch->strategy) {
193:   case SCOTCH_STRATQUALITY:     *strategy = MP_PTSCOTCH_QUALITY; break;
194:   case SCOTCH_STRATSPEED:       *strategy = MP_PTSCOTCH_SPEED; break;
195:   case SCOTCH_STRATBALANCE:     *strategy = MP_PTSCOTCH_BALANCE; break;
196:   case SCOTCH_STRATSAFETY:      *strategy = MP_PTSCOTCH_SAFETY; break;
197:   case SCOTCH_STRATSCALABILITY: *strategy = MP_PTSCOTCH_SCALABILITY; break;
198:   }
199:   return(0);
200: }

204: PetscErrorCode MatPartitioningView_PTScotch(MatPartitioning part, PetscViewer viewer)
205: {
206:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
207:   PetscErrorCode           ierr;
208:   PetscBool                isascii;
209:   const char               *str=0;

212:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
213:   if (isascii) {
214:     switch (scotch->strategy) {
215:     case SCOTCH_STRATQUALITY:     str = "Prioritize quality over speed"; break;
216:     case SCOTCH_STRATSPEED:       str = "Prioritize speed over quality"; break;
217:     case SCOTCH_STRATBALANCE:     str = "Enforce load balance"; break;
218:     case SCOTCH_STRATSAFETY:      str = "Avoid methods that may fail"; break;
219:     case SCOTCH_STRATSCALABILITY: str = "Favor scalability as much as possible"; break;
220:     }
221:     PetscViewerASCIIPrintf(viewer,"  Strategy=%s\n",str);
222:     PetscViewerASCIIPrintf(viewer,"  Load imbalance ratio=%g\n",scotch->imbalance);
223:   }
224:   return(0);
225: }

229: PetscErrorCode MatPartitioningSetFromOptions_PTScotch(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
230: {
231:   PetscErrorCode           ierr;
232:   PetscBool                flag;
233:   PetscReal                r;
234:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
235:   MPPTScotchStrategyType   strat;

238:   MatPartitioningPTScotchGetStrategy(part,&strat);
239:   PetscOptionsHead(PetscOptionsObject,"PTScotch partitioning options");
240:   PetscOptionsEnum("-mat_partitioning_ptscotch_strategy","Strategy","MatPartitioningPTScotchSetStrategy",MPPTScotchStrategyTypes,(PetscEnum)strat,(PetscEnum*)&strat,&flag);
241:   if (flag) { MatPartitioningPTScotchSetStrategy(part,strat); }
242:   PetscOptionsReal("-mat_partitioning_ptscotch_imbalance","Load imbalance ratio","MatPartitioningPTScotchSetImbalance",scotch->imbalance,&r,&flag);
243:   if (flag) { MatPartitioningPTScotchSetImbalance(part,r); }
244:   PetscOptionsTail();
245:   return(0);
246: }

250: PetscErrorCode MatPartitioningApply_PTScotch(MatPartitioning part,IS *partitioning)
251: {
252:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
253:   PetscErrorCode           ierr;
254:   PetscMPIInt              rank;
255:   Mat                      mat  = part->adj;
256:   Mat_MPIAdj               *adj = (Mat_MPIAdj*)mat->data;
257:   PetscBool                flg;
258:   PetscInt                 i,j,wgtflag=0,bs=1,nold;
259:   PetscReal                *vwgttab,deltval;
260:   SCOTCH_Num               *locals,*velotab,*veloloctab,*edloloctab,vertlocnbr,edgelocnbr,nparts=part->n;
261:   SCOTCH_Arch              archdat;
262:   SCOTCH_Dgraph            grafdat;
263:   SCOTCH_Dmapping          mappdat;
264:   SCOTCH_Strat             stradat;

267:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
268:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
269:   if (!flg) {
270:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
271:        resulting partition results need to be stretched to match the original matrix */
272:     nold = mat->rmap->n;
273:     MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat);
274:     bs   = nold/mat->rmap->n;
275:     adj  = (Mat_MPIAdj*)mat->data;
276:   }

278:   PetscMalloc1(mat->rmap->n+1,&locals);
279:   PetscMalloc1(nparts,&vwgttab);
280:   PetscMalloc1(nparts,&velotab);
281:   for (j=0; j<nparts; j++) {
282:     if (part->part_weights) vwgttab[j] = part->part_weights[j]*nparts;
283:     else vwgttab[j] = 1.0;
284:   }
285:   for (i=0; i<nparts; i++) {
286:     deltval = PetscAbsReal(vwgttab[i]-PetscFloorReal(vwgttab[i]+0.5));
287:     if (deltval>0.01) {
288:       for (j=0; j<nparts; j++) vwgttab[j] /= deltval;
289:     }
290:   }
291:   for (i=0; i<nparts; i++) velotab[i] = (SCOTCH_Num)(vwgttab[i] + 0.5);
292:   PetscFree(vwgttab);

294:   SCOTCH_dgraphInit(&grafdat,PetscObjectComm((PetscObject)part));

296:   vertlocnbr = mat->rmap->range[rank+1] - mat->rmap->range[rank];
297:   edgelocnbr = adj->i[vertlocnbr];
298:   veloloctab = (!part->vertex_weights && !(wgtflag & 2)) ? part->vertex_weights : NULL;
299:   edloloctab = (!adj->values && !(wgtflag & 1)) ? adj->values : NULL;

301:   SCOTCH_dgraphBuild(&grafdat,0,vertlocnbr,vertlocnbr,adj->i,adj->i+1,veloloctab,
302:                             NULL,edgelocnbr,edgelocnbr,adj->j,NULL,edloloctab);

304: #if defined(PETSC_USE_DEBUG)
305:   SCOTCH_dgraphCheck(&grafdat);
306: #endif

308:   SCOTCH_archInit(&archdat);
309:   SCOTCH_stratInit(&stradat);
310:   SCOTCH_stratDgraphMapBuild(&stradat,scotch->strategy,nparts,nparts,scotch->imbalance);

312:   SCOTCH_archCmpltw(&archdat,nparts,velotab);
313:   SCOTCH_dgraphMapInit(&grafdat,&mappdat,&archdat,locals);
314:   SCOTCH_dgraphMapCompute(&grafdat,&mappdat,&stradat);

316:   SCOTCH_dgraphMapExit (&grafdat,&mappdat);
317:   SCOTCH_archExit(&archdat);
318:   SCOTCH_stratExit(&stradat);
319:   SCOTCH_dgraphExit(&grafdat);
320:   PetscFree(velotab);

322:   if (bs > 1) {
323:     PetscInt *newlocals;
324:     PetscMalloc1(bs*mat->rmap->n,&newlocals);
325:     for (i=0;i<mat->rmap->n;i++) {
326:       for (j=0;j<bs;j++) {
327:         newlocals[bs*i+j] = locals[i];
328:       }
329:     }
330:     PetscFree(locals);
331:     ISCreateGeneral(PetscObjectComm((PetscObject)part),bs*mat->rmap->n,newlocals,PETSC_OWN_POINTER,partitioning);
332:   } else {
333:     ISCreateGeneral(PetscObjectComm((PetscObject)part),mat->rmap->n,locals,PETSC_OWN_POINTER,partitioning);
334:   }

336:   if (!flg) {
337:     MatDestroy(&mat);
338:   }
339:   return(0);
340: }

344: PetscErrorCode MatPartitioningDestroy_PTScotch(MatPartitioning part)
345: {
346:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
347:   PetscErrorCode           ierr;

350:   PetscFree(scotch);
351:   /* clear composed functions */
352:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetImbalance_C",NULL);
353:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetImbalance_C",NULL);
354:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetStrategy_C",NULL);
355:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetStrategy_C",NULL);
356:   return(0);
357: }

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

362:    Level: beginner

364:    Notes: See http://www.labri.fr/perso/pelegrin/scotch/

366: .keywords: Partitioning, create, context

368: .seealso: MatPartitioningSetType(), MatPartitioningType
369: M*/

373: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning part)
374: {
375:   PetscErrorCode           ierr;
376:   MatPartitioning_PTScotch *scotch;

379:   PetscNewLog(part,&scotch);
380:   part->data = (void*)scotch;

382:   scotch->imbalance = 0.01;
383:   scotch->strategy  = SCOTCH_STRATQUALITY;

385:   part->ops->apply          = MatPartitioningApply_PTScotch;
386:   part->ops->view           = MatPartitioningView_PTScotch;
387:   part->ops->setfromoptions = MatPartitioningSetFromOptions_PTScotch;
388:   part->ops->destroy        = MatPartitioningDestroy_PTScotch;

390:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetImbalance_C",MatPartitioningPTScotchSetImbalance_PTScotch);
391:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetImbalance_C",MatPartitioningPTScotchGetImbalance_PTScotch);
392:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetStrategy_C",MatPartitioningPTScotchSetStrategy_PTScotch);
393:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetStrategy_C",MatPartitioningPTScotchGetStrategy_PTScotch);
394:   return(0);
395: }