Actual source code: partition.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <petsc/private/matimpl.h>

  4: /* Logging support */
  5: PetscClassId MAT_PARTITIONING_CLASSID;

  7: /*
  8:    Simplest partitioning, keeps the current partitioning.
  9: */
 10: static PetscErrorCode MatPartitioningApply_Current(MatPartitioning part,IS *partitioning)
 11: {
 13:   PetscInt       m;
 14:   PetscMPIInt    rank,size;

 17:   MPI_Comm_size(PetscObjectComm((PetscObject)part),&size);
 18:   if (part->n != size) {
 19:     const char *prefix;
 20:     PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
 21:     SETERRQ1(PetscObjectComm((PetscObject)part),PETSC_ERR_SUP,"This is the DEFAULT NO-OP partitioner, it currently only supports one domain per processor\nuse -%smat_partitioning_type parmetis or chaco or ptscotch for more than one subdomain per processor",prefix ? prefix : "");
 22:   }
 23:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);

 25:   MatGetLocalSize(part->adj,&m,NULL);
 26:   ISCreateStride(PetscObjectComm((PetscObject)part),m,rank,0,partitioning);
 27:   return(0);
 28: }

 30: /*
 31:    partition an index to rebalance the computation
 32: */
 33: static PetscErrorCode MatPartitioningApply_Average(MatPartitioning part,IS *partitioning)
 34: {
 36:   PetscInt       m,M,nparts,*indices,r,d,*parts,i,start,end,loc;

 39:   MatGetSize(part->adj,&M,NULL);
 40:   MatGetLocalSize(part->adj,&m,NULL);
 41:   nparts = part->n;
 42:   PetscCalloc1(nparts,&parts);
 43:   d = M/nparts;
 44:   for (i=0; i<nparts; i++) parts[i] = d;
 45:   r = M%nparts;
 46:   for (i=0; i<r; i++) parts[i] += 1;
 47:   for (i=1; i<nparts; i++) parts[i] += parts[i-1];
 48:   PetscCalloc1(m,&indices);
 49:   MatGetOwnershipRange(part->adj,&start,&end);
 50:   for (i=start; i<end; i++) {
 51:     PetscFindInt(i,nparts,parts,&loc);
 52:     if (loc<0) loc = -(loc+1);
 53:     else loc = loc+1;
 54:     indices[i-start] = loc;
 55:   }
 56:   PetscFree(parts);
 57:   ISCreateGeneral(PetscObjectComm((PetscObject)part),m,indices,PETSC_OWN_POINTER,partitioning);
 58:   return(0);
 59: }

 61: static PetscErrorCode MatPartitioningApply_Square(MatPartitioning part,IS *partitioning)
 62: {
 64:   PetscInt       cell,n,N,p,rstart,rend,*color;
 65:   PetscMPIInt    size;

 68:   MPI_Comm_size(PetscObjectComm((PetscObject)part),&size);
 69:   if (part->n != size) SETERRQ(PetscObjectComm((PetscObject)part),PETSC_ERR_SUP,"Currently only supports one domain per processor");
 70:   p = (PetscInt)PetscSqrtReal((PetscReal)part->n);
 71:   if (p*p != part->n) SETERRQ(PetscObjectComm((PetscObject)part),PETSC_ERR_SUP,"Square partitioning requires \"perfect square\" number of domains");

 73:   MatGetSize(part->adj,&N,NULL);
 74:   n    = (PetscInt)PetscSqrtReal((PetscReal)N);
 75:   if (n*n != N) SETERRQ(PetscObjectComm((PetscObject)part),PETSC_ERR_SUP,"Square partitioning requires square domain");
 76:   if (n%p != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Square partitioning requires p to divide n");
 77:   MatGetOwnershipRange(part->adj,&rstart,&rend);
 78:   PetscMalloc1(rend-rstart,&color);
 79:   /* for (int cell=rstart; cell<rend; cell++) { color[cell-rstart] = ((cell%n) < (n/2)) + 2 * ((cell/n) < (n/2)); } */
 80:   for (cell=rstart; cell<rend; cell++) {
 81:     color[cell-rstart] = ((cell%n) / (n/p)) + p * ((cell/n) / (n/p));
 82:   }
 83:   ISCreateGeneral(PetscObjectComm((PetscObject)part),rend-rstart,color,PETSC_OWN_POINTER,partitioning);
 84:   return(0);
 85: }

 87: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Current(MatPartitioning part)
 88: {
 90:   part->ops->apply   = MatPartitioningApply_Current;
 91:   part->ops->view    = 0;
 92:   part->ops->destroy = 0;
 93:   return(0);
 94: }

 96: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Average(MatPartitioning part)
 97: {
 99:   part->ops->apply   = MatPartitioningApply_Average;
100:   part->ops->view    = 0;
101:   part->ops->destroy = 0;
102:   return(0);
103: }

105: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Square(MatPartitioning part)
106: {
108:   part->ops->apply   = MatPartitioningApply_Square;
109:   part->ops->view    = 0;
110:   part->ops->destroy = 0;
111:   return(0);
112: }


115: /* gets as input the "sizes" array computed by ParMetis_*_NodeND and returns
116:        seps[  0 :         2*p) : the start and end node of each subdomain
117:        seps[2*p : 2*p+2*(p-1)) : the start and end node of each separator
118:      levels[  0 :         p-1) : level in the tree for each separator (-1 root, -2 and -3 first level and so on)
119:    The arrays must be large enough
120: */
121: PETSC_INTERN PetscErrorCode MatPartitioningSizesToSep_Private(PetscInt p, PetscInt sizes[], PetscInt seps[], PetscInt level[])
122: {
123:   PetscInt       l2p,i,pTree,pStartTree;

127:   l2p = PetscLog2Real(p);
128:   if (l2p - (PetscInt)PetscLog2Real(p)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"%D is not a power of 2",p);
129:   if (!p) return(0);
130:   PetscMemzero(seps,(2*p-2)*sizeof(PetscInt));
131:   PetscMemzero(level,(p-1)*sizeof(PetscInt));
132:   seps[2*p-2] = sizes[2*p-2];
133:   pTree = p;
134:   pStartTree = 0;
135:   while (pTree != 1) {
136:     for (i = pStartTree; i < pStartTree + pTree; i++) {
137:       seps[i] += sizes[i];
138:       seps[pStartTree + pTree + (i-pStartTree)/2] += seps[i];
139:     }
140:     pStartTree += pTree;
141:     pTree = pTree/2;
142:   }
143:   seps[2*p-2] -= sizes[2*p-2];

145:   pStartTree = 2*p-2;
146:   pTree      = 1;
147:   while (pStartTree > 0) {
148:     for (i = pStartTree; i < pStartTree + pTree; i++) {
149:       PetscInt k = 2*i - (pStartTree +2*pTree);
150:       PetscInt n = seps[k+1];

152:       seps[k+1]  = seps[i]   - sizes[k+1];
153:       seps[k]    = seps[k+1] + sizes[k+1] - n - sizes[k];
154:       level[i-p] = -pTree - i + pStartTree;
155:     }
156:     pTree *= 2;
157:     pStartTree -= pTree;
158:   }
159:   /* I know there should be a formula */
160:   PetscSortIntWithArrayPair(p-1,seps+p,sizes+p,level);
161:   for (i=2*p-2;i>=0;i--) { seps[2*i] = seps[i]; seps[2*i+1] = seps[i] + sizes[i] - 1; }
162:   return(0);
163: }

165: /* ===========================================================================================*/

167: PetscFunctionList MatPartitioningList              = 0;
168: PetscBool         MatPartitioningRegisterAllCalled = PETSC_FALSE;


171: /*@C
172:    MatPartitioningRegister - Adds a new sparse matrix partitioning to the  matrix package.

174:    Not Collective

176:    Input Parameters:
177: +  sname - name of partitioning (for example MATPARTITIONINGCURRENT) or parmetis
178: -  function - function pointer that creates the partitioning type

180:    Level: developer

182:    Sample usage:
183: .vb
184:    MatPartitioningRegister("my_part",MyPartCreate);
185: .ve

187:    Then, your partitioner can be chosen with the procedural interface via
188: $     MatPartitioningSetType(part,"my_part")
189:    or at runtime via the option
190: $     -mat_partitioning_type my_part

192: .keywords: matrix, partitioning, register

194: .seealso: MatPartitioningRegisterDestroy(), MatPartitioningRegisterAll()
195: @*/
196: PetscErrorCode  MatPartitioningRegister(const char sname[],PetscErrorCode (*function)(MatPartitioning))
197: {

201:   MatInitializePackage();
202:   PetscFunctionListAdd(&MatPartitioningList,sname,function);
203:   return(0);
204: }

206: /*@C
207:    MatPartitioningGetType - Gets the Partitioning method type and name (as a string)
208:         from the partitioning context.

210:    Not collective

212:    Input Parameter:
213: .  partitioning - the partitioning context

215:    Output Parameter:
216: .  type - partitioner type

218:    Level: intermediate

220:    Not Collective

222: .keywords: Partitioning, get, method, name, type
223: @*/
224: PetscErrorCode  MatPartitioningGetType(MatPartitioning partitioning,MatPartitioningType *type)
225: {
229:   *type = ((PetscObject)partitioning)->type_name;
230:   return(0);
231: }

233: /*@C
234:    MatPartitioningSetNParts - Set how many partitions need to be created;
235:         by default this is one per processor. Certain partitioning schemes may
236:         in fact only support that option.

238:    Not collective

240:    Input Parameter:
241: .  partitioning - the partitioning context
242: .  n - the number of partitions

244:    Level: intermediate

246:    Not Collective

248: .keywords: Partitioning, set

250: .seealso: MatPartitioningCreate(), MatPartitioningApply()
251: @*/
252: PetscErrorCode  MatPartitioningSetNParts(MatPartitioning part,PetscInt n)
253: {
255:   part->n = n;
256:   return(0);
257: }

259: /*@
260:    MatPartitioningApplyND - Gets a nested dissection partitioning for a matrix.

262:    Collective on Mat

264:    Input Parameters:
265: .  matp - the matrix partitioning object

267:    Output Parameters:
268: .   partitioning - the partitioning. For each local node, a positive value indicates the processor
269:                    number the node has been assigned to. Negative x values indicate the separator level -(x+1).

271:    Level: beginner

273:    The user can define additional partitionings; see MatPartitioningRegister().

275: .keywords: matrix, get, partitioning

277: .seealso:  MatPartitioningRegister(), MatPartitioningCreate(),
278:            MatPartitioningDestroy(), MatPartitioningSetAdjacency(), ISPartitioningToNumbering(),
279:            ISPartitioningCount()
280: @*/
281: PetscErrorCode  MatPartitioningApplyND(MatPartitioning matp,IS *partitioning)
282: {

288:   if (!matp->adj->assembled) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
289:   if (matp->adj->factortype) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
290:   if (!matp->ops->applynd) SETERRQ1(PetscObjectComm((PetscObject)matp),PETSC_ERR_SUP,"Nested dissection not provided by MatPartitioningType %s",((PetscObject)matp)->type_name);
291:   PetscLogEventBegin(MAT_PartitioningND,matp,0,0,0);
292:   (*matp->ops->applynd)(matp,partitioning);
293:   PetscLogEventEnd(MAT_PartitioningND,matp,0,0,0);

295:   MatPartitioningViewFromOptions(matp,NULL,"-mat_partitioning_view");
296:   ISViewFromOptions(*partitioning,NULL,"-mat_partitioning_view");
297:   return(0);
298: }

300: /*@
301:    MatPartitioningApply - Gets a partitioning for a matrix.

303:    Collective on Mat

305:    Input Parameters:
306: .  matp - the matrix partitioning object

308:    Output Parameters:
309: .   partitioning - the partitioning. For each local node this tells the processor
310:                    number that that node is assigned to.

312:    Options Database Keys:
313:    To specify the partitioning through the options database, use one of
314:    the following
315: $    -mat_partitioning_type parmetis, -mat_partitioning current
316:    To see the partitioning result
317: $    -mat_partitioning_view

319:    Level: beginner

321:    The user can define additional partitionings; see MatPartitioningRegister().

323: .keywords: matrix, get, partitioning

325: .seealso:  MatPartitioningRegister(), MatPartitioningCreate(),
326:            MatPartitioningDestroy(), MatPartitioningSetAdjacency(), ISPartitioningToNumbering(),
327:            ISPartitioningCount()
328: @*/
329: PetscErrorCode  MatPartitioningApply(MatPartitioning matp,IS *partitioning)
330: {
332:   PetscBool      viewbalance,improve;

337:   if (!matp->adj->assembled) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
338:   if (matp->adj->factortype) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
339:   if (!matp->ops->apply) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Must set type with MatPartitioningSetFromOptions() or MatPartitioningSetType()");
340:   PetscLogEventBegin(MAT_Partitioning,matp,0,0,0);
341:   (*matp->ops->apply)(matp,partitioning);
342:   PetscLogEventEnd(MAT_Partitioning,matp,0,0,0);

344:   MatPartitioningViewFromOptions(matp,NULL,"-mat_partitioning_view");
345:   ISViewFromOptions(*partitioning,NULL,"-mat_partitioning_view");

347:   PetscObjectOptionsBegin((PetscObject)matp);
348:   viewbalance = PETSC_FALSE;
349:   PetscOptionsBool("-mat_partitioning_view_imbalance","Display imbalance information of a partition",NULL,PETSC_FALSE,&viewbalance,NULL);
350:   improve = PETSC_FALSE;
351:   PetscOptionsBool("-mat_partitioning_improve","Improve the quality of a partition",NULL,PETSC_FALSE,&improve,NULL);
352:   PetscOptionsEnd();

354:   if (improve) {
355:     MatPartitioningImprove(matp,partitioning);
356:   }

358:   if (viewbalance) {
359:     MatPartitioningViewImbalance(matp,*partitioning);
360:   }
361:   return(0);
362: }


365: /*@
366:    MatPartitioningImprove - Improves the quality of a given partition.

368:    Collective on Mat

370:    Input Parameters:
371: .  matp - the matrix partitioning object
372: .  partitioning - the partitioning. For each local node this tells the processor
373:                    number that that node is assigned to.

375:    Output Parameters:
376: .   partitioning - the partitioning. For each local node this tells the processor
377:                    number that that node is assigned to.

379:    Options Database Keys:
380:    To improve the quality of the partition
381: $    -mat_partitioning_improve

383:    Level: beginner


386: .keywords: matrix, improve, partitioning

388: .seealso:  MatPartitioningApply(), MatPartitioningCreate(),
389:            MatPartitioningDestroy(), MatPartitioningSetAdjacency(), ISPartitioningToNumbering(),
390:            ISPartitioningCount()
391: @*/
392: PetscErrorCode  MatPartitioningImprove(MatPartitioning matp,IS *partitioning)
393: {

399:   if (!matp->adj->assembled) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix");
400:   if (matp->adj->factortype) SETERRQ(PetscObjectComm((PetscObject)matp),PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix");
401:   PetscLogEventBegin(MAT_Partitioning,matp,0,0,0);
402:   if (matp->ops->improve) {
403:     (*matp->ops->improve)(matp,partitioning);
404:   }
405:   PetscLogEventEnd(MAT_Partitioning,matp,0,0,0);
406:   return(0);
407: }

409: /*@
410:    MatPartitioningViewImbalance - Display partitioning imbalance information.

412:    Collective on MatPartitioning

414:    Input Parameters:
415: .  matp - the matrix partitioning object
416: .  partitioning - the partitioning. For each local node this tells the processor
417:                    number that that node is assigned to.

419:    Options Database Keys:
420:    To see the partitioning imbalance information
421: $    -mat_partitioning_view_balance

423:    Level: beginner

425: .keywords: matrix, imbalance, partitioning

427: .seealso:  MatPartitioningApply(), MatPartitioningView()
428: @*/
429: PetscErrorCode  MatPartitioningViewImbalance(MatPartitioning matp, IS partitioning)
430: {
431:   PetscErrorCode  ierr;
432:   PetscInt        nparts,*subdomainsizes,*subdomainsizes_tmp,nlocal,i,maxsub,minsub,avgsub;
433:   const PetscInt  *indices;
434:   PetscViewer     viewer;

439:   nparts = matp->n;
440:   PetscCalloc2(nparts,&subdomainsizes,nparts,&subdomainsizes_tmp);
441:   ISGetLocalSize(partitioning,&nlocal);
442:   ISGetIndices(partitioning,&indices);
443:   for (i=0;i<nlocal;i++) {
444:     subdomainsizes_tmp[indices[i]] += matp->vertex_weights? matp->vertex_weights[i]:1;
445:   }
446:   MPI_Allreduce(subdomainsizes_tmp,subdomainsizes,nparts,MPIU_INT,MPI_SUM, PetscObjectComm((PetscObject)matp));
447:   ISRestoreIndices(partitioning,&indices);
448:   minsub = PETSC_MAX_INT, maxsub = PETSC_MIN_INT, avgsub=0;
449:   for (i=0; i<nparts; i++) {
450:     minsub = PetscMin(minsub,subdomainsizes[i]);
451:     maxsub = PetscMax(maxsub,subdomainsizes[i]);
452:     avgsub += subdomainsizes[i];
453:   }
454:   avgsub /=nparts;
455:   PetscFree2(subdomainsizes,subdomainsizes_tmp);
456:   PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)matp),&viewer);
457:   MatPartitioningView(matp,viewer);
458:   PetscViewerASCIIPrintf(viewer,"Partitioning Imbalance Info: Max %D, Min %D, Avg %D, R %g\n",maxsub, minsub, avgsub, (double)(maxsub/(PetscReal)minsub));
459:   return(0);
460: }

462: /*@
463:    MatPartitioningSetAdjacency - Sets the adjacency graph (matrix) of the thing to be
464:       partitioned.

466:    Collective on MatPartitioning and Mat

468:    Input Parameters:
469: +  part - the partitioning context
470: -  adj - the adjacency matrix

472:    Level: beginner

474: .keywords: Partitioning, adjacency

476: .seealso: MatPartitioningCreate()
477: @*/
478: PetscErrorCode  MatPartitioningSetAdjacency(MatPartitioning part,Mat adj)
479: {
483:   part->adj = adj;
484:   return(0);
485: }

487: /*@
488:    MatPartitioningDestroy - Destroys the partitioning context.

490:    Collective on Partitioning

492:    Input Parameters:
493: .  part - the partitioning context

495:    Level: beginner

497: .keywords: Partitioning, destroy, context

499: .seealso: MatPartitioningCreate()
500: @*/
501: PetscErrorCode  MatPartitioningDestroy(MatPartitioning *part)
502: {

506:   if (!*part) return(0);
508:   if (--((PetscObject)(*part))->refct > 0) {*part = 0; return(0);}

510:   if ((*part)->ops->destroy) {
511:     (*(*part)->ops->destroy)((*part));
512:   }
513:   PetscFree((*part)->vertex_weights);
514:   PetscFree((*part)->part_weights);
515:   PetscHeaderDestroy(part);
516:   return(0);
517: }

519: /*@C
520:    MatPartitioningSetVertexWeights - Sets the weights for vertices for a partitioning.

522:    Logically Collective on Partitioning

524:    Input Parameters:
525: +  part - the partitioning context
526: -  weights - the weights, on each process this array must have the same size as the number of local rows

528:    Level: beginner

530:    Notes:
531:       The array weights is freed by PETSc so the user should not free the array. In C/C++
532:    the array must be obtained with a call to PetscMalloc(), not malloc().

534: .keywords: Partitioning, destroy, context

536: .seealso: MatPartitioningCreate(), MatPartitioningSetType(), MatPartitioningSetPartitionWeights()
537: @*/
538: PetscErrorCode  MatPartitioningSetVertexWeights(MatPartitioning part,const PetscInt weights[])
539: {


545:   PetscFree(part->vertex_weights);

547:   part->vertex_weights = (PetscInt*)weights;
548:   return(0);
549: }

551: /*@C
552:    MatPartitioningSetPartitionWeights - Sets the weights for each partition.

554:    Logically Collective on Partitioning

556:    Input Parameters:
557: +  part - the partitioning context
558: -  weights - An array of size nparts that is used to specify the fraction of
559:              vertex weight that should be distributed to each sub-domain for
560:              the balance constraint. If all of the sub-domains are to be of
561:              the same size, then each of the nparts elements should be set
562:              to a value of 1/nparts. Note that the sum of all of the weights
563:              should be one.

565:    Level: beginner

567:    Notes:
568:       The array weights is freed by PETSc so the user should not free the array. In C/C++
569:    the array must be obtained with a call to PetscMalloc(), not malloc().

571: .keywords: Partitioning, destroy, context

573: .seealso: MatPartitioningCreate(), MatPartitioningSetType(), MatPartitioningSetVertexWeights()
574: @*/
575: PetscErrorCode  MatPartitioningSetPartitionWeights(MatPartitioning part,const PetscReal weights[])
576: {


582:   PetscFree(part->part_weights);

584:   part->part_weights = (PetscReal*)weights;
585:   return(0);
586: }

588: /*@
589:    MatPartitioningCreate - Creates a partitioning context.

591:    Collective on MPI_Comm

593:    Input Parameter:
594: .   comm - MPI communicator

596:    Output Parameter:
597: .  newp - location to put the context

599:    Level: beginner

601: .keywords: Partitioning, create, context

603: .seealso: MatPartitioningSetType(), MatPartitioningApply(), MatPartitioningDestroy(),
604:           MatPartitioningSetAdjacency()

606: @*/
607: PetscErrorCode  MatPartitioningCreate(MPI_Comm comm,MatPartitioning *newp)
608: {
609:   MatPartitioning part;
610:   PetscErrorCode  ierr;
611:   PetscMPIInt     size;

614:   *newp = 0;

616:   MatInitializePackage();
617:   PetscHeaderCreate(part,MAT_PARTITIONING_CLASSID,"MatPartitioning","Matrix/graph partitioning","MatOrderings",comm,MatPartitioningDestroy,MatPartitioningView);
618:   part->vertex_weights = NULL;
619:   part->part_weights   = NULL;

621:   MPI_Comm_size(comm,&size);
622:   part->n = (PetscInt)size;

624:   *newp = part;
625:   return(0);
626: }

628: /*@C
629:    MatPartitioningView - Prints the partitioning data structure.

631:    Collective on MatPartitioning

633:    Input Parameters:
634: .  part - the partitioning context
635: .  viewer - optional visualization context

637:    Level: intermediate

639:    Note:
640:    The available visualization contexts include
641: +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
642: -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
643:          output where only the first processor opens
644:          the file.  All other processors send their
645:          data to the first processor to print.

647:    The user can open alternative visualization contexts with
648: .     PetscViewerASCIIOpen() - output to a specified file

650: .keywords: Partitioning, view

652: .seealso: PetscViewerASCIIOpen()
653: @*/
654: PetscErrorCode  MatPartitioningView(MatPartitioning part,PetscViewer viewer)
655: {
657:   PetscBool      iascii;

661:   if (!viewer) {
662:     PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)part),&viewer);
663:   }

667:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
668:   if (iascii) {
669:     PetscObjectPrintClassNamePrefixType((PetscObject)part,viewer);
670:     if (part->vertex_weights) {
671:       PetscViewerASCIIPrintf(viewer,"  Using vertex weights\n");
672:     }
673:   }
674:   if (part->ops->view) {
675:     PetscViewerASCIIPushTab(viewer);
676:     (*part->ops->view)(part,viewer);
677:     PetscViewerASCIIPopTab(viewer);
678:   }
679:   return(0);
680: }

682: /*@C
683:    MatPartitioningSetType - Sets the type of partitioner to use

685:    Collective on MatPartitioning

687:    Input Parameter:
688: .  part - the partitioning context.
689: .  type - a known method

691:    Options Database Command:
692: $  -mat_partitioning_type  <type>
693: $      Use -help for a list of available methods
694: $      (for instance, parmetis)

696:    Level: intermediate

698: .keywords: partitioning, set, method, type

700: .seealso: MatPartitioningCreate(), MatPartitioningApply(), MatPartitioningType

702: @*/
703: PetscErrorCode  MatPartitioningSetType(MatPartitioning part,MatPartitioningType type)
704: {
705:   PetscErrorCode ierr,(*r)(MatPartitioning);
706:   PetscBool      match;


712:   PetscObjectTypeCompare((PetscObject)part,type,&match);
713:   if (match) return(0);

715:   if (part->ops->destroy) {
716:     (*part->ops->destroy)(part);
717:     part->ops->destroy = NULL;
718:   }
719:   part->setupcalled = 0;
720:   part->data        = 0;
721:   PetscMemzero(part->ops,sizeof(struct _MatPartitioningOps));

723:   PetscFunctionListFind(MatPartitioningList,type,&r);
724:   if (!r) SETERRQ1(PetscObjectComm((PetscObject)part),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown partitioning type %s",type);

726:   (*r)(part);

728:   PetscFree(((PetscObject)part)->type_name);
729:   PetscStrallocpy(type,&((PetscObject)part)->type_name);
730:   return(0);
731: }

733: /*@
734:    MatPartitioningSetFromOptions - Sets various partitioning options from the
735:         options database.

737:    Collective on MatPartitioning

739:    Input Parameter:
740: .  part - the partitioning context.

742:    Options Database Command:
743: $  -mat_partitioning_type  <type>
744: $      Use -help for a list of available methods
745: $      (for instance, parmetis)
746: $  -mat_partitioning_nparts - number of subgraphs


749:    Notes:
750:     If the partitioner has not been set by the user it uses one of the installed partitioner such as ParMetis. If there are
751:    no installed partitioners it uses current which means no repartioning.

753:    Level: beginner

755: .keywords: partitioning, set, method, type
756: @*/
757: PetscErrorCode  MatPartitioningSetFromOptions(MatPartitioning part)
758: {
760:   PetscBool      flag;
761:   char           type[256];
762:   const char     *def;

765:   PetscObjectOptionsBegin((PetscObject)part);
766:   if (!((PetscObject)part)->type_name) {
767: #if defined(PETSC_HAVE_PARMETIS)
768:     def = MATPARTITIONINGPARMETIS;
769: #elif defined(PETSC_HAVE_CHACO)
770:     def = MATPARTITIONINGCHACO;
771: #elif defined(PETSC_HAVE_PARTY)
772:     def = MATPARTITIONINGPARTY;
773: #elif defined(PETSC_HAVE_PTSCOTCH)
774:     def = MATPARTITIONINGPTSCOTCH;
775: #else
776:     def = MATPARTITIONINGCURRENT;
777: #endif
778:   } else {
779:     def = ((PetscObject)part)->type_name;
780:   }
781:   PetscOptionsFList("-mat_partitioning_type","Type of partitioner","MatPartitioningSetType",MatPartitioningList,def,type,256,&flag);
782:   if (flag) {
783:     MatPartitioningSetType(part,type);
784:   }

786:   PetscOptionsInt("-mat_partitioning_nparts","number of fine parts",NULL,part->n,& part->n,&flag);

788:   /*
789:     Set the type if it was never set.
790:   */
791:   if (!((PetscObject)part)->type_name) {
792:     MatPartitioningSetType(part,def);
793:   }

795:   if (part->ops->setfromoptions) {
796:     (*part->ops->setfromoptions)(PetscOptionsObject,part);
797:   }
798:   PetscOptionsEnd();
799:   return(0);
800: }