Actual source code: scotch.c

petsc-3.9.4 2018-09-11
Report Typos and Errors

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

  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;

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

 17:    Collective on MatPartitioning

 19:    Input Parameters:
 20: +  part - the partitioning context
 21: -  imb  - the load imbalance ratio

 23:    Options Database:
 24: .  -mat_partitioning_ptscotch_imbalance <imb>

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

 29:    Level: advanced

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

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

 44: PetscErrorCode MatPartitioningPTScotchSetImbalance_PTScotch(MatPartitioning part,PetscReal imb)
 45: {
 46:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

 49:   if (imb==PETSC_DEFAULT) scotch->imbalance = 0.01;
 50:   else {
 51:     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]");
 52:     scotch->imbalance = (double)imb;
 53:   }
 54:   return(0);
 55: }

 57: /*@
 58:    MatPartitioningPTScotchGetImbalance - Gets the value of the load imbalance
 59:    ratio used during strategy selection.

 61:    Not Collective

 63:    Input Parameter:
 64: .  part - the partitioning context

 66:    Output Parameter:
 67: .  imb  - the load imbalance ratio

 69:    Level: advanced

 71: .seealso: MatPartitioningPTScotchSetImbalance()
 72: @*/
 73: PetscErrorCode MatPartitioningPTScotchGetImbalance(MatPartitioning part,PetscReal *imb)
 74: {

 80:   PetscUseMethod(part,"MatPartitioningPTScotchGetImbalance_C",(MatPartitioning,PetscReal*),(part,imb));
 81:   return(0);
 82: }

 84: PetscErrorCode MatPartitioningPTScotchGetImbalance_PTScotch(MatPartitioning part,PetscReal *imb)
 85: {
 86:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

 89:   *imb = scotch->imbalance;
 90:   return(0);
 91: }

 93: /*@
 94:    MatPartitioningPTScotchSetStrategy - Sets the strategy to be used in PTScotch.

 96:    Collective on MatPartitioning

 98:    Input Parameters:
 99: +  part - the partitioning context
100: -  strategy - the strategy, one of
101: .vb
102:      MP_PTSCOTCH_DEFAULT     - Default behavior
103:      MP_PTSCOTCH_QUALITY     - Prioritize quality over speed
104:      MP_PTSCOTCH_SPEED       - Prioritize speed over quality
105:      MP_PTSCOTCH_BALANCE     - Enforce load balance
106:      MP_PTSCOTCH_SAFETY      - Avoid methods that may fail
107:      MP_PTSCOTCH_SCALABILITY - Favor scalability as much as possible
108: .ve

110:    Options Database:
111: .  -mat_partitioning_ptscotch_strategy [quality,speed,balance,safety,scalability] - strategy

113:    Level: advanced

115:    Notes:
116:    The default is MP_SCOTCH_QUALITY. See the PTScotch documentation for more information.

118: .seealso: MatPartitioningPTScotchSetImbalance(), MatPartitioningPTScotchGetStrategy()
119: @*/
120: PetscErrorCode MatPartitioningPTScotchSetStrategy(MatPartitioning part,MPPTScotchStrategyType strategy)
121: {

127:   PetscTryMethod(part,"MatPartitioningPTScotchSetStrategy_C",(MatPartitioning,MPPTScotchStrategyType),(part,strategy));
128:   return(0);
129: }

131: PetscErrorCode MatPartitioningPTScotchSetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType strategy)
132: {
133:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

136:   switch (strategy) {
137:   case MP_PTSCOTCH_QUALITY:     scotch->strategy = SCOTCH_STRATQUALITY; break;
138:   case MP_PTSCOTCH_SPEED:       scotch->strategy = SCOTCH_STRATSPEED; break;
139:   case MP_PTSCOTCH_BALANCE:     scotch->strategy = SCOTCH_STRATBALANCE; break;
140:   case MP_PTSCOTCH_SAFETY:      scotch->strategy = SCOTCH_STRATSAFETY; break;
141:   case MP_PTSCOTCH_SCALABILITY: scotch->strategy = SCOTCH_STRATSCALABILITY; break;
142:   default:                      scotch->strategy = SCOTCH_STRATDEFAULT; break;
143:   }
144:   return(0);
145: }

147: /*@
148:    MatPartitioningPTScotchGetStrategy - Gets the strategy used in PTScotch.

150:    Not Collective

152:    Input Parameter:
153: .  part - the partitioning context

155:    Output Parameter:
156: .  strategy - the strategy

158:    Level: advanced

160: .seealso: MatPartitioningPTScotchSetStrategy()
161: @*/
162: PetscErrorCode MatPartitioningPTScotchGetStrategy(MatPartitioning part,MPPTScotchStrategyType *strategy)
163: {

169:   PetscUseMethod(part,"MatPartitioningPTScotchGetStrategy_C",(MatPartitioning,MPPTScotchStrategyType*),(part,strategy));
170:   return(0);
171: }

173: PetscErrorCode MatPartitioningPTScotchGetStrategy_PTScotch(MatPartitioning part,MPPTScotchStrategyType *strategy)
174: {
175:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;

178:   switch (scotch->strategy) {
179:   case SCOTCH_STRATQUALITY:     *strategy = MP_PTSCOTCH_QUALITY; break;
180:   case SCOTCH_STRATSPEED:       *strategy = MP_PTSCOTCH_SPEED; break;
181:   case SCOTCH_STRATBALANCE:     *strategy = MP_PTSCOTCH_BALANCE; break;
182:   case SCOTCH_STRATSAFETY:      *strategy = MP_PTSCOTCH_SAFETY; break;
183:   case SCOTCH_STRATSCALABILITY: *strategy = MP_PTSCOTCH_SCALABILITY; break;
184:   default:                      *strategy = MP_PTSCOTCH_DEFAULT; break;
185:   }
186:   return(0);
187: }

189: PetscErrorCode MatPartitioningView_PTScotch(MatPartitioning part, PetscViewer viewer)
190: {
191:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
192:   PetscErrorCode           ierr;
193:   PetscBool                isascii;
194:   const char               *str=0;

197:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
198:   if (isascii) {
199:     switch (scotch->strategy) {
200:     case SCOTCH_STRATQUALITY:     str = "Prioritize quality over speed"; break;
201:     case SCOTCH_STRATSPEED:       str = "Prioritize speed over quality"; break;
202:     case SCOTCH_STRATBALANCE:     str = "Enforce load balance"; break;
203:     case SCOTCH_STRATSAFETY:      str = "Avoid methods that may fail"; break;
204:     case SCOTCH_STRATSCALABILITY: str = "Favor scalability as much as possible"; break;
205:     default:                      str = "Default behavior"; break;
206:     }
207:     PetscViewerASCIIPrintf(viewer,"  Strategy=%s\n",str);
208:     PetscViewerASCIIPrintf(viewer,"  Load imbalance ratio=%g\n",scotch->imbalance);
209:   }
210:   return(0);
211: }

213: PetscErrorCode MatPartitioningSetFromOptions_PTScotch(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
214: {
215:   PetscErrorCode           ierr;
216:   PetscBool                flag;
217:   PetscReal                r;
218:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
219:   MPPTScotchStrategyType   strat;

222:   MatPartitioningPTScotchGetStrategy(part,&strat);
223:   PetscOptionsHead(PetscOptionsObject,"PTScotch partitioning options");
224:   PetscOptionsEnum("-mat_partitioning_ptscotch_strategy","Strategy","MatPartitioningPTScotchSetStrategy",MPPTScotchStrategyTypes,(PetscEnum)strat,(PetscEnum*)&strat,&flag);
225:   if (flag) { MatPartitioningPTScotchSetStrategy(part,strat); }
226:   PetscOptionsReal("-mat_partitioning_ptscotch_imbalance","Load imbalance ratio","MatPartitioningPTScotchSetImbalance",scotch->imbalance,&r,&flag);
227:   if (flag) { MatPartitioningPTScotchSetImbalance(part,r); }
228:   PetscOptionsTail();
229:   return(0);
230: }

232: PetscErrorCode MatPartitioningApply_PTScotch(MatPartitioning part,IS *partitioning)
233: {
234:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
235:   PetscErrorCode           ierr;
236:   PetscMPIInt              rank;
237:   Mat                      mat  = part->adj;
238:   Mat_MPIAdj               *adj = (Mat_MPIAdj*)mat->data;
239:   PetscBool                flg,distributed;
240:   PetscBool                proc_weight_flg;
241:   PetscInt                 i,j,p,bs=1,nold;
242:   PetscReal                *vwgttab,deltval;
243:   SCOTCH_Num               *locals,*velotab,*veloloctab,*edloloctab,vertlocnbr,edgelocnbr,nparts=part->n;

246:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
247:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
248:   if (!flg) {
249:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
250:        resulting partition results need to be stretched to match the original matrix */
251:     nold = mat->rmap->n;
252:     MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&mat);
253:     if (mat->rmap->n > 0) bs = nold/mat->rmap->n;
254:     adj  = (Mat_MPIAdj*)mat->data;
255:   }

257:   proc_weight_flg = PETSC_TRUE;
258:   PetscOptionsGetBool(NULL, NULL, "-mat_partitioning_ptscotch_proc_weight", &proc_weight_flg, NULL);

260:   PetscMalloc1(mat->rmap->n+1,&locals);

262:   velotab = NULL;
263:   if (proc_weight_flg)
264:   {
265:   PetscMalloc1(nparts,&vwgttab);
266:   PetscMalloc1(nparts,&velotab);
267:   for (j=0; j<nparts; j++) {
268:     if (part->part_weights) vwgttab[j] = part->part_weights[j]*nparts;
269:     else vwgttab[j] = 1.0;
270:   }
271:   for (i=0; i<nparts; i++) {
272:     deltval = PetscAbsReal(vwgttab[i]-PetscFloorReal(vwgttab[i]+0.5));
273:     if (deltval>0.01) {
274:       for (j=0; j<nparts; j++) vwgttab[j] /= deltval;
275:     }
276:   }
277:   for (i=0; i<nparts; i++) velotab[i] = (SCOTCH_Num)(vwgttab[i] + 0.5);
278:   PetscFree(vwgttab);
279:   }

281:   vertlocnbr = mat->rmap->range[rank+1] - mat->rmap->range[rank];
282:   edgelocnbr = adj->i[vertlocnbr];
283:   veloloctab = part->vertex_weights;
284:   edloloctab = adj->values;

286:   /* detect whether all vertices are located at the same process in original graph */
287:   for (p = 0; !mat->rmap->range[p+1] && p < nparts; ++p);
288:   distributed = (mat->rmap->range[p+1] == mat->rmap->N) ? PETSC_FALSE : PETSC_TRUE;

290:   if (distributed) {
291:     SCOTCH_Arch              archdat;
292:     SCOTCH_Dgraph            grafdat;
293:     SCOTCH_Dmapping          mappdat;
294:     SCOTCH_Strat             stradat;

296:     SCOTCH_dgraphInit(&grafdat,PetscObjectComm((PetscObject)part));
297:     SCOTCH_dgraphBuild(&grafdat,0,vertlocnbr,vertlocnbr,adj->i,adj->i+1,veloloctab,
298:                               NULL,edgelocnbr,edgelocnbr,adj->j,NULL,edloloctab);

300: #if defined(PETSC_USE_DEBUG)
301:     SCOTCH_dgraphCheck(&grafdat);
302: #endif

304:     SCOTCH_archInit(&archdat);
305:     SCOTCH_stratInit(&stradat);
306:     SCOTCH_stratDgraphMapBuild(&stradat,scotch->strategy,nparts,nparts,scotch->imbalance);

308:     if (velotab) {
309:       SCOTCH_archCmpltw(&archdat,nparts,velotab);
310:     } else {
311:       SCOTCH_archCmplt( &archdat,nparts);
312:     }
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);

321:   } else if (rank == p) {
322:     SCOTCH_Graph   grafdat;
323:     SCOTCH_Strat   stradat;

325:     SCOTCH_graphInit(&grafdat);
326:     SCOTCH_graphBuild(&grafdat,0,vertlocnbr,adj->i,adj->i+1,veloloctab,NULL,edgelocnbr,adj->j,edloloctab);

328: #if defined(PETSC_USE_DEBUG)
329:     SCOTCH_graphCheck(&grafdat);
330: #endif

332:     SCOTCH_stratInit(&stradat);
333:     SCOTCH_stratGraphMapBuild(&stradat,scotch->strategy,nparts,scotch->imbalance);

335:     SCOTCH_graphPart(&grafdat,nparts,&stradat,locals);

337:     SCOTCH_stratExit(&stradat);
338:     SCOTCH_graphExit(&grafdat);
339:   }

341:   PetscFree(velotab);

343:   if (bs > 1) {
344:     PetscInt *newlocals;
345:     PetscMalloc1(bs*mat->rmap->n,&newlocals);
346:     for (i=0;i<mat->rmap->n;i++) {
347:       for (j=0;j<bs;j++) {
348:         newlocals[bs*i+j] = locals[i];
349:       }
350:     }
351:     PetscFree(locals);
352:     ISCreateGeneral(PetscObjectComm((PetscObject)part),bs*mat->rmap->n,newlocals,PETSC_OWN_POINTER,partitioning);
353:   } else {
354:     ISCreateGeneral(PetscObjectComm((PetscObject)part),mat->rmap->n,locals,PETSC_OWN_POINTER,partitioning);
355:   }

357:   if (!flg) {
358:     MatDestroy(&mat);
359:   }
360:   return(0);
361: }

363: PetscErrorCode MatPartitioningDestroy_PTScotch(MatPartitioning part)
364: {
365:   MatPartitioning_PTScotch *scotch = (MatPartitioning_PTScotch*)part->data;
366:   PetscErrorCode           ierr;

369:   PetscFree(scotch);
370:   /* clear composed functions */
371:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetImbalance_C",NULL);
372:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetImbalance_C",NULL);
373:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetStrategy_C",NULL);
374:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetStrategy_C",NULL);
375:   return(0);
376: }

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

381:    Level: beginner

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

385: .keywords: Partitioning, create, context

387: .seealso: MatPartitioningSetType(), MatPartitioningType
388: M*/

390: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_PTScotch(MatPartitioning part)
391: {
392:   PetscErrorCode           ierr;
393:   MatPartitioning_PTScotch *scotch;

396:   PetscNewLog(part,&scotch);
397:   part->data = (void*)scotch;

399:   scotch->imbalance = 0.01;
400:   scotch->strategy  = SCOTCH_STRATDEFAULT;

402:   part->ops->apply          = MatPartitioningApply_PTScotch;
403:   part->ops->view           = MatPartitioningView_PTScotch;
404:   part->ops->setfromoptions = MatPartitioningSetFromOptions_PTScotch;
405:   part->ops->destroy        = MatPartitioningDestroy_PTScotch;

407:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetImbalance_C",MatPartitioningPTScotchSetImbalance_PTScotch);
408:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetImbalance_C",MatPartitioningPTScotchGetImbalance_PTScotch);
409:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchSetStrategy_C",MatPartitioningPTScotchSetStrategy_PTScotch);
410:   PetscObjectComposeFunction((PetscObject)part,"MatPartitioningPTScotchGetStrategy_C",MatPartitioningPTScotchGetStrategy_PTScotch);
411:   return(0);
412: }