Actual source code: pmetis.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

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

  4: /*
  5:    Currently using ParMetis-4.0.2
  6: */

  8: #include <parmetis.h>

 10: /*
 11:       The first 5 elements of this structure are the input control array to Metis
 12: */
 13: typedef struct {
 14:   PetscInt  cuts;         /* number of cuts made (output) */
 15:   PetscInt  foldfactor;
 16:   PetscInt  parallel;     /* use parallel partitioner for coarse problem */
 17:   PetscInt  indexing;     /* 0 indicates C indexing, 1 Fortran */
 18:   PetscInt  printout;     /* indicates if one wishes Metis to print info */
 19:   PetscBool repartition;
 20: } MatPartitioning_Parmetis;

 22: #define CHKERRQPARMETIS(n,func)                                             \
 23:   if (n == METIS_ERROR_INPUT) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to wrong inputs and/or options for %s",func); \
 24:   else if (n == METIS_ERROR_MEMORY) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS error due to insufficient memory in %s",func); \
 25:   else if (n == METIS_ERROR) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"ParMETIS general error in %s",func); \

 27: #define PetscStackCallParmetis(func,args) do {PetscStackPush(#func);int status = func args;PetscStackPop; CHKERRQPARMETIS(status,#func);} while (0)

 29: static PetscErrorCode MatPartitioningApply_Parmetis_Private(MatPartitioning part, PetscBool useND, PetscBool isImprove, IS *partitioning)
 30: {
 31:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*)part->data;
 32:   PetscErrorCode           ierr;
 33:   PetscInt                 *locals = NULL;
 34:   Mat                      mat     = part->adj,amat,pmat;
 35:   PetscBool                flg;
 36:   PetscInt                 bs = 1;

 41:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
 42:   if (flg) {
 43:     amat = mat;
 44:     PetscObjectReference((PetscObject)amat);
 45:   } else {
 46:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
 47:        resulting partition results need to be stretched to match the original matrix */
 48:     MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&amat);
 49:     if (amat->rmap->n > 0) bs = mat->rmap->n/amat->rmap->n;
 50:   }
 51:   MatMPIAdjCreateNonemptySubcommMat(amat,&pmat);
 52:   MPI_Barrier(PetscObjectComm((PetscObject)part));

 54:   if (pmat) {
 55:     MPI_Comm   pcomm,comm;
 56:     Mat_MPIAdj *adj     = (Mat_MPIAdj*)pmat->data;
 57:     PetscInt   *vtxdist = pmat->rmap->range;
 58:     PetscInt   *xadj    = adj->i;
 59:     PetscInt   *adjncy  = adj->j;
 60:     PetscInt   *NDorder = NULL;
 61:     PetscInt   itmp     = 0,wgtflag=0, numflag=0, ncon=1, nparts=part->n, options[24], i, j;
 62:     real_t     *tpwgts,*ubvec,itr=0.1;

 64:     PetscObjectGetComm((PetscObject)pmat,&pcomm);
 65: #if defined(PETSC_USE_DEBUG)
 66:     /* check that matrix has no diagonal entries */
 67:     {
 68:       PetscInt rstart;
 69:       MatGetOwnershipRange(pmat,&rstart,NULL);
 70:       for (i=0; i<pmat->rmap->n; i++) {
 71:         for (j=xadj[i]; j<xadj[i+1]; j++) {
 72:           if (adjncy[j] == i+rstart) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Row %D has diagonal entry; Parmetis forbids diagonal entry",i+rstart);
 73:         }
 74:       }
 75:     }
 76: #endif

 78:     PetscMalloc1(pmat->rmap->n,&locals);

 80:     if (isImprove) {
 81:       PetscInt       i;
 82:       const PetscInt *part_indices;
 84:       ISGetIndices(*partitioning,&part_indices);
 85:       for (i=0; i<pmat->rmap->n; i++) locals[i] = part_indices[i*bs];
 86:       ISRestoreIndices(*partitioning,&part_indices);
 87:       ISDestroy(partitioning);
 88:     }

 90:     if (adj->values && part->use_edge_weights && !part->vertex_weights) wgtflag = 1;
 91:     if (part->vertex_weights && !adj->values) wgtflag = 2;
 92:     if (part->vertex_weights && adj->values && part->use_edge_weights) wgtflag = 3;

 94:     if (PetscLogPrintInfo) {itmp = pmetis->printout; pmetis->printout = 127;}
 95:     PetscMalloc1(ncon*nparts,&tpwgts);
 96:     for (i=0; i<ncon; i++) {
 97:       for (j=0; j<nparts; j++) {
 98:         if (part->part_weights) {
 99:           tpwgts[i*nparts+j] = part->part_weights[i*nparts+j];
100:         } else {
101:           tpwgts[i*nparts+j] = 1./nparts;
102:         }
103:       }
104:     }
105:     PetscMalloc1(ncon,&ubvec);
106:     for (i=0; i<ncon; i++) ubvec[i] = 1.05;
107:     /* This sets the defaults */
108:     options[0] = 0;
109:     for (i=1; i<24; i++) options[i] = -1;
110:     /* Duplicate the communicator to be sure that ParMETIS attribute caching does not interfere with PETSc. */
111:     MPI_Comm_dup(pcomm,&comm);
112:     if (useND) {
113:       PetscInt    *sizes, *seps, log2size, subd, *level;
114:       PetscMPIInt size;
115:       idx_t       mtype = PARMETIS_MTYPE_GLOBAL, rtype = PARMETIS_SRTYPE_2PHASE, p_nseps = 1, s_nseps = 1;
116:       real_t      ubfrac = 1.05;

118:       MPI_Comm_size(comm,&size);
119:       PetscMalloc1(pmat->rmap->n,&NDorder);
120:       PetscMalloc3(2*size,&sizes,4*size,&seps,size,&level);
121:       PetscStackCallParmetis(ParMETIS_V32_NodeND,((idx_t*)vtxdist,(idx_t*)xadj,(idx_t*)adjncy,(idx_t*)part->vertex_weights,(idx_t*)&numflag,&mtype,&rtype,&p_nseps,&s_nseps,&ubfrac,NULL/* seed */,NULL/* dbglvl */,(idx_t*)NDorder,(idx_t*)(sizes),&comm));
122:       log2size = PetscLog2Real(size);
123:       subd = PetscPowInt(2,log2size);
124:       MatPartitioningSizesToSep_Private(subd,sizes,seps,level);
125:       for (i=0;i<pmat->rmap->n;i++) {
126:         PetscInt loc;

128:         PetscFindInt(NDorder[i],2*subd,seps,&loc);
129:         if (loc < 0) {
130:           loc = -(loc+1);
131:           if (loc%2) { /* part of subdomain */
132:             locals[i] = loc/2;
133:           } else {
134:             PetscFindInt(NDorder[i],2*(subd-1),seps+2*subd,&loc);
135:             loc = loc < 0 ? -(loc+1)/2 : loc/2;
136:             locals[i] = level[loc];
137:           }
138:         } else locals[i] = loc/2;
139:       }
140:       PetscFree3(sizes,seps,level);
141:     } else {
142:       if (pmetis->repartition) {
143:         PetscStackCallParmetis(ParMETIS_V3_AdaptiveRepart,((idx_t*)vtxdist,(idx_t*)xadj,(idx_t*)adjncy,(idx_t*)part->vertex_weights,(idx_t*)part->vertex_weights,(idx_t*)adj->values,(idx_t*)&wgtflag,(idx_t*)&numflag,(idx_t*)&ncon,(idx_t*)&nparts,tpwgts,ubvec,&itr,(idx_t*)options,(idx_t*)&pmetis->cuts,(idx_t*)locals,&comm));
144:       } else if (isImprove) {
145:         PetscStackCallParmetis(ParMETIS_V3_RefineKway,((idx_t*)vtxdist,(idx_t*)xadj,(idx_t*)adjncy,(idx_t*)part->vertex_weights,(idx_t*)adj->values,(idx_t*)&wgtflag,(idx_t*)&numflag,(idx_t*)&ncon,(idx_t*)&nparts,tpwgts,ubvec,(idx_t*)options,(idx_t*)&pmetis->cuts,(idx_t*)locals,&comm));
146:       } else {
147:         PetscStackCallParmetis(ParMETIS_V3_PartKway,((idx_t*)vtxdist,(idx_t*)xadj,(idx_t*)adjncy,(idx_t*)part->vertex_weights,(idx_t*)adj->values,(idx_t*)&wgtflag,(idx_t*)&numflag,(idx_t*)&ncon,(idx_t*)&nparts,tpwgts,ubvec,(idx_t*)options,(idx_t*)&pmetis->cuts,(idx_t*)locals,&comm));
148:       }
149:     }
150:     MPI_Comm_free(&comm);

152:     PetscFree(tpwgts);
153:     PetscFree(ubvec);
154:     if (PetscLogPrintInfo) pmetis->printout = itmp;

156:     if (bs > 1) {
157:       PetscInt i,j,*newlocals;
158:       PetscMalloc1(bs*pmat->rmap->n,&newlocals);
159:       for (i=0; i<pmat->rmap->n; i++) {
160:         for (j=0; j<bs; j++) {
161:           newlocals[bs*i + j] = locals[i];
162:         }
163:       }
164:       PetscFree(locals);
165:       ISCreateGeneral(PetscObjectComm((PetscObject)part),bs*pmat->rmap->n,newlocals,PETSC_OWN_POINTER,partitioning);
166:     } else {
167:       ISCreateGeneral(PetscObjectComm((PetscObject)part),pmat->rmap->n,locals,PETSC_OWN_POINTER,partitioning);
168:     }
169:     if (useND) {
170:       IS ndis;

172:       if (bs > 1) {
173:         ISCreateBlock(PetscObjectComm((PetscObject)part),bs,pmat->rmap->n,NDorder,PETSC_OWN_POINTER,&ndis);
174:       } else {
175:         ISCreateGeneral(PetscObjectComm((PetscObject)part),pmat->rmap->n,NDorder,PETSC_OWN_POINTER,&ndis);
176:       }
177:       ISSetPermutation(ndis);
178:       PetscObjectCompose((PetscObject)(*partitioning),"_petsc_matpartitioning_ndorder",(PetscObject)ndis);
179:       ISDestroy(&ndis);
180:     }
181:   } else {
182:     ISCreateGeneral(PetscObjectComm((PetscObject)part),0,NULL,PETSC_COPY_VALUES,partitioning);
183:     if (useND) {
184:       IS ndis;

186:       if (bs > 1) {
187:         ISCreateBlock(PetscObjectComm((PetscObject)part),bs,0,NULL,PETSC_COPY_VALUES,&ndis);
188:       } else {
189:         ISCreateGeneral(PetscObjectComm((PetscObject)part),0,NULL,PETSC_COPY_VALUES,&ndis);
190:       }
191:       ISSetPermutation(ndis);
192:       PetscObjectCompose((PetscObject)(*partitioning),"_petsc_matpartitioning_ndorder",(PetscObject)ndis);
193:       ISDestroy(&ndis);
194:     }
195:   }
196:   MatDestroy(&pmat);
197:   MatDestroy(&amat);
198:   return(0);
199: }

201: /*
202:    Uses the ParMETIS parallel matrix partitioner to compute a nested dissection ordering of the matrix in parallel
203: */
204: static PetscErrorCode MatPartitioningApplyND_Parmetis(MatPartitioning part, IS *partitioning)
205: {

209:   MatPartitioningApply_Parmetis_Private(part, PETSC_TRUE, PETSC_FALSE, partitioning);
210:   return(0);
211: }

213: /*
214:    Uses the ParMETIS parallel matrix partitioner to partition the matrix in parallel
215: */
216: static PetscErrorCode MatPartitioningApply_Parmetis(MatPartitioning part, IS *partitioning)
217: {

221:   MatPartitioningApply_Parmetis_Private(part, PETSC_FALSE, PETSC_FALSE, partitioning);
222:   return(0);
223: }

225: /*
226:    Uses the ParMETIS to improve the quality  of a partition
227: */
228: static PetscErrorCode MatPartitioningImprove_Parmetis(MatPartitioning part, IS *partitioning)
229: {

233:   MatPartitioningApply_Parmetis_Private(part, PETSC_FALSE, PETSC_TRUE, partitioning);
234:   return(0);
235: }

237: PetscErrorCode MatPartitioningView_Parmetis(MatPartitioning part,PetscViewer viewer)
238: {
239:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*)part->data;
240:   PetscErrorCode           ierr;
241:   PetscMPIInt              rank;
242:   PetscBool                iascii;

245:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
246:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
247:   if (iascii) {
248:     if (pmetis->parallel == 2) {
249:       PetscViewerASCIIPrintf(viewer,"  Using parallel coarse grid partitioner\n");
250:     } else {
251:       PetscViewerASCIIPrintf(viewer,"  Using sequential coarse grid partitioner\n");
252:     }
253:     PetscViewerASCIIPrintf(viewer,"  Using %D fold factor\n",pmetis->foldfactor);
254:     PetscViewerASCIIPushSynchronized(viewer);
255:     PetscViewerASCIISynchronizedPrintf(viewer,"  [%d]Number of cuts found %D\n",rank,pmetis->cuts);
256:     PetscViewerFlush(viewer);
257:     PetscViewerASCIIPopSynchronized(viewer);
258:   }
259:   return(0);
260: }

262: /*@
263:      MatPartitioningParmetisSetCoarseSequential - Use the sequential code to
264:          do the partitioning of the coarse grid.

266:   Logically Collective on MatPartitioning

268:   Input Parameter:
269: .  part - the partitioning context

271:    Level: advanced

273: @*/
274: PetscErrorCode  MatPartitioningParmetisSetCoarseSequential(MatPartitioning part)
275: {
276:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*)part->data;

279:   pmetis->parallel = 1;
280:   return(0);
281: }

283: /*@
284:      MatPartitioningParmetisSetRepartition - Repartition
285:      current mesh to rebalance computation.

287:   Logically Collective on MatPartitioning

289:   Input Parameter:
290: .  part - the partitioning context

292:    Level: advanced

294: @*/
295: PetscErrorCode  MatPartitioningParmetisSetRepartition(MatPartitioning part)
296: {
297:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*)part->data;

300:   pmetis->repartition = PETSC_TRUE;
301:   return(0);
302: }

304: /*@
305:   MatPartitioningParmetisGetEdgeCut - Returns the number of edge cuts in the vertex partition.

307:   Input Parameter:
308: . part - the partitioning context

310:   Output Parameter:
311: . cut - the edge cut

313:    Level: advanced

315: @*/
316: PetscErrorCode  MatPartitioningParmetisGetEdgeCut(MatPartitioning part, PetscInt *cut)
317: {
318:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*) part->data;

321:   *cut = pmetis->cuts;
322:   return(0);
323: }

325: PetscErrorCode MatPartitioningSetFromOptions_Parmetis(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
326: {
328:   PetscBool      flag = PETSC_FALSE;

331:   PetscOptionsHead(PetscOptionsObject,"Set ParMeTiS partitioning options");
332:   PetscOptionsBool("-mat_partitioning_parmetis_coarse_sequential","Use sequential coarse partitioner","MatPartitioningParmetisSetCoarseSequential",flag,&flag,NULL);
333:   if (flag) {
334:     MatPartitioningParmetisSetCoarseSequential(part);
335:   }
336:   PetscOptionsBool("-mat_partitioning_parmetis_repartition","","MatPartitioningParmetisSetRepartition",flag,&flag,NULL);
337:   if(flag){
338:      MatPartitioningParmetisSetRepartition(part);
339:   }
340:   PetscOptionsTail();
341:   return(0);
342: }


345: PetscErrorCode MatPartitioningDestroy_Parmetis(MatPartitioning part)
346: {
347:   MatPartitioning_Parmetis *pmetis = (MatPartitioning_Parmetis*)part->data;
348:   PetscErrorCode           ierr;

351:   PetscFree(pmetis);
352:   return(0);
353: }


356: /*MC
357:    MATPARTITIONINGPARMETIS - Creates a partitioning context via the external package PARMETIS.

359:    Collective

361:    Input Parameter:
362: .  part - the partitioning context

364:    Options Database Keys:
365: .  -mat_partitioning_parmetis_coarse_sequential - use sequential PARMETIS coarse partitioner

367:    Level: beginner

369:    Notes:
370:     See https://www-users.cs.umn.edu/~karypis/metis/

372: .seealso: MatPartitioningSetType(), MatPartitioningType

374: M*/

376: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Parmetis(MatPartitioning part)
377: {
378:   PetscErrorCode           ierr;
379:   MatPartitioning_Parmetis *pmetis;

382:   PetscNewLog(part,&pmetis);
383:   part->data = (void*)pmetis;

385:   pmetis->cuts       = 0;   /* output variable */
386:   pmetis->foldfactor = 150; /*folding factor */
387:   pmetis->parallel   = 2;   /* use parallel partitioner for coarse grid */
388:   pmetis->indexing   = 0;   /* index numbering starts from 0 */
389:   pmetis->printout   = 0;   /* print no output while running */
390:   pmetis->repartition      = PETSC_FALSE;

392:   part->ops->apply          = MatPartitioningApply_Parmetis;
393:   part->ops->applynd        = MatPartitioningApplyND_Parmetis;
394:   part->ops->improve        = MatPartitioningImprove_Parmetis;
395:   part->ops->view           = MatPartitioningView_Parmetis;
396:   part->ops->destroy        = MatPartitioningDestroy_Parmetis;
397:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Parmetis;
398:   return(0);
399: }

401: /*@
402:  MatMeshToVertexGraph -   This routine does not exist because ParMETIS does not provide the functionality.  Uses the ParMETIS package to
403:                        convert a Mat that represents a mesh to a Mat the represents the graph of the coupling
404:                        between vertices of the cells and is suitable for partitioning with the MatPartitioning object. Use this to partition
405:                        vertices of a mesh. More likely you should use MatMeshToCellGraph()

407:    Collective on Mat

409:    Input Parameter:
410: +     mesh - the graph that represents the mesh
411: -     ncommonnodes - mesh elements that share this number of common nodes are considered neighbors, use 2 for triangles and
412:                      quadrilaterials, 3 for tetrahedrals and 4 for hexahedrals

414:    Output Parameter:
415: .     dual - the dual graph

417:    Notes:
418:      Currently requires ParMetis to be installed and uses ParMETIS_V3_Mesh2Dual()

420:      The columns of each row of the Mat mesh are the global vertex numbers of the vertices of that rows cell. The number of rows in mesh is
421:      number of cells, the number of columns is the number of vertices.

423:    Level: advanced

425: .seealso: MatMeshToCellGraph(), MatCreateMPIAdj(), MatPartitioningCreate()

427: @*/
428: PetscErrorCode MatMeshToVertexGraph(Mat mesh,PetscInt ncommonnodes,Mat *dual)
429: {
431:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"ParMETIS does not provide this functionality");
432:   return(0);
433: }

435: /*@
436:      MatMeshToCellGraph -   Uses the ParMETIS package to convert a Mat that represents a mesh to a Mat the represents the graph of the coupling
437:                        between cells (the "dual" graph) and is suitable for partitioning with the MatPartitioning object. Use this to partition
438:                        cells of a mesh.

440:    Collective on Mat

442:    Input Parameter:
443: +     mesh - the graph that represents the mesh
444: -     ncommonnodes - mesh elements that share this number of common nodes are considered neighbors, use 2 for triangles and
445:                      quadrilaterials, 3 for tetrahedrals and 4 for hexahedrals

447:    Output Parameter:
448: .     dual - the dual graph

450:    Notes:
451:      Currently requires ParMetis to be installed and uses ParMETIS_V3_Mesh2Dual()

453: $     Each row of the mesh object represents a single cell in the mesh. For triangles it has 3 entries, quadrilaterials 4 entries,
454: $         tetrahedrals 4 entries and hexahedrals 8 entries. You can mix triangles and quadrilaterals in the same mesh, but cannot
455: $         mix  tetrahedrals and hexahedrals
456: $     The columns of each row of the Mat mesh are the global vertex numbers of the vertices of that row's cell.
457: $     The number of rows in mesh is number of cells, the number of columns is the number of vertices.


460:    Level: advanced

462: .seealso: MatMeshToVertexGraph(), MatCreateMPIAdj(), MatPartitioningCreate()


465: @*/
466: PetscErrorCode MatMeshToCellGraph(Mat mesh,PetscInt ncommonnodes,Mat *dual)
467: {
469:   PetscInt       *newxadj,*newadjncy;
470:   PetscInt       numflag=0;
471:   Mat_MPIAdj     *adj   = (Mat_MPIAdj*)mesh->data,*newadj;
472:   PetscBool      flg;
473:   MPI_Comm       comm;

476:   PetscObjectTypeCompare((PetscObject)mesh,MATMPIADJ,&flg);
477:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Must use MPIAdj matrix type");

479:   PetscObjectGetComm((PetscObject)mesh,&comm);
480:   PetscStackCallParmetis(ParMETIS_V3_Mesh2Dual,((idx_t*)mesh->rmap->range,(idx_t*)adj->i,(idx_t*)adj->j,(idx_t*)&numflag,(idx_t*)&ncommonnodes,(idx_t**)&newxadj,(idx_t**)&newadjncy,&comm));
481:   MatCreateMPIAdj(PetscObjectComm((PetscObject)mesh),mesh->rmap->n,mesh->rmap->N,newxadj,newadjncy,NULL,dual);
482:   newadj = (Mat_MPIAdj*)(*dual)->data;

484:   newadj->freeaijwithfree = PETSC_TRUE; /* signal the matrix should be freed with system free since space was allocated by ParMETIS */
485:   return(0);
486: }