Actual source code: hierarchical.c


  2: #include <../src/mat/impls/adj/mpi/mpiadj.h>
  3: #include <petscsf.h>
  4: #include <petsc/private/matimpl.h>

  6: /*
  7:   It is a hierarchical partitioning. The partitioner has two goals:
  8:   (1) Most of current partitioners fail at a large scale. The hierarchical partitioning
  9:   strategy is trying to produce large number of subdomains when number of processor cores is large.
 10:   (2) PCGASM needs one 'big' subdomain across multi-cores. The partitioner provides two
 11:   consistent partitions, coarse parts and fine parts. A coarse part is a 'big' subdomain consisting
 12:   of several small subdomains.
 13: */

 15: PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning,IS,PetscInt,PetscInt,IS*);
 16: PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat,IS,IS,IS*,Mat*,ISLocalToGlobalMapping*);
 17: PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat,IS,ISLocalToGlobalMapping,IS*);

 19: typedef struct {
 20:   char*                fineparttype; /* partitioner on fine level */
 21:   char*                coarseparttype; /* partitioner on coarse level */
 22:   PetscInt             nfineparts; /* number of fine parts on each coarse subdomain */
 23:   PetscInt             ncoarseparts; /* number of coarse parts */
 24:   IS                   coarseparts; /* partitioning on coarse level */
 25:   IS                   fineparts; /* partitioning on fine level */
 26:   MatPartitioning      coarseMatPart; /* MatPartititioning on coarse level (first level) */
 27:   MatPartitioning      fineMatPart; /* MatPartitioning on fine level (second level) */
 28:   MatPartitioning      improver; /* Improve the quality of a partition */
 29: } MatPartitioning_Hierarchical;

 31: /*
 32:    Uses a hierarchical partitioning strategy to partition the matrix in parallel.
 33:    Use this interface to make the partitioner consistent with others
 34: */
 35: static PetscErrorCode MatPartitioningApply_Hierarchical(MatPartitioning part,IS *partitioning)
 36: {
 37:   MatPartitioning_Hierarchical *hpart  = (MatPartitioning_Hierarchical*)part->data;
 38:   const PetscInt               *fineparts_indices, *coarseparts_indices;
 39:   PetscInt                     *fineparts_indices_tmp;
 40:   PetscInt                     *parts_indices,i,j,mat_localsize, *offsets;
 41:   Mat                           mat    = part->adj,adj,sadj;
 42:   PetscReal                    *part_weights;
 43:   PetscBool                     flg;
 44:   PetscInt                      bs     = 1;
 45:   PetscInt                     *coarse_vertex_weights = NULL;
 46:   PetscMPIInt                   size,rank;
 47:   MPI_Comm                      comm,scomm;
 48:   IS                            destination,fineparts_temp, vweights, svweights;
 49:   PetscInt                      nsvwegihts,*fp_vweights;
 50:   const PetscInt                *svweights_indices;
 51:   ISLocalToGlobalMapping        mapping;
 52:   const char                    *prefix;
 53:   PetscBool                     use_edge_weights;

 55:   PetscObjectGetComm((PetscObject)part,&comm);
 56:   MPI_Comm_size(comm,&size);
 57:   MPI_Comm_rank(comm,&rank);
 58:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
 59:   if (flg) {
 60:     adj = mat;
 61:     PetscObjectReference((PetscObject)adj);
 62:   }else {
 63:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
 64:        resulting partition results need to be stretched to match the original matrix */
 65:    MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
 66:    if (adj->rmap->n > 0) bs = mat->rmap->n/adj->rmap->n;
 67:   }
 68:   /* local size of mat */
 69:   mat_localsize = adj->rmap->n;
 70:   /* check parameters */
 71:   /* how many small subdomains we want from a given 'big' suddomain */

 75:   /* Partitioning the domain into one single subdomain is a trivial case, and we should just return  */
 76:   if (part->n==1) {
 77:     PetscCalloc1(bs*adj->rmap->n,&parts_indices);
 78:     ISCreateGeneral(comm,bs*adj->rmap->n,parts_indices,PETSC_OWN_POINTER,partitioning);
 79:     hpart->ncoarseparts = 1;
 80:     hpart->nfineparts = 1;
 81:     PetscStrallocpy("NONE",&hpart->coarseparttype);
 82:     PetscStrallocpy("NONE",&hpart->fineparttype);
 83:     MatDestroy(&adj);
 84:     return 0;
 85:   }

 87:   if (part->n) {
 88:     hpart->ncoarseparts = part->n/hpart->nfineparts;

 90:     if (part->n%hpart->nfineparts != 0) hpart->ncoarseparts++;
 91:   }else{
 92:     part->n = hpart->ncoarseparts*hpart->nfineparts;
 93:   }

 95:   PetscMalloc1(hpart->ncoarseparts+1, &offsets);
 96:   PetscMalloc1(hpart->ncoarseparts, &part_weights);

 98:   offsets[0] = 0;
 99:   if (part->n%hpart->nfineparts != 0) offsets[1] = part->n%hpart->nfineparts;
100:   else offsets[1] = hpart->nfineparts;

102:   part_weights[0] = ((PetscReal)offsets[1])/part->n;

104:   for (i=2; i<=hpart->ncoarseparts; i++) {
105:     offsets[i] = hpart->nfineparts;
106:     part_weights[i-1] = ((PetscReal)offsets[i])/part->n;
107:   }

109:   offsets[0] = 0;
110:   for (i=1;i<=hpart->ncoarseparts; i++)
111:     offsets[i] += offsets[i-1];

113:   /* If these exists a mat partitioner, we should delete it */
114:   MatPartitioningDestroy(&hpart->coarseMatPart);
115:   MatPartitioningCreate(comm,&hpart->coarseMatPart);
116:   PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
117:   PetscObjectSetOptionsPrefix((PetscObject)hpart->coarseMatPart,prefix);
118:   PetscObjectAppendOptionsPrefix((PetscObject)hpart->coarseMatPart,"hierarch_coarse_");
119:     /* if did not set partitioning type yet, use parmetis by default */
120:   if (!hpart->coarseparttype) {
121: #if defined(PETSC_HAVE_PARMETIS)
122:     MatPartitioningSetType(hpart->coarseMatPart,MATPARTITIONINGPARMETIS);
123:     PetscStrallocpy(MATPARTITIONINGPARMETIS,&hpart->coarseparttype);
124: #elif defined(PETSC_HAVE_PTSCOTCH)
125:     MatPartitioningSetType(hpart->coarseMatPart,MATPARTITIONINGPTSCOTCH);
126:     PetscStrallocpy(MATPARTITIONINGPTSCOTCH,&hpart->coarseparttype);
127: #else
128:     SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
129: #endif
130:   } else {
131:     MatPartitioningSetType(hpart->coarseMatPart,hpart->coarseparttype);
132:   }
133:   MatPartitioningSetAdjacency(hpart->coarseMatPart,adj);
134:   MatPartitioningSetNParts(hpart->coarseMatPart, hpart->ncoarseparts);
135:   /* copy over vertex weights */
136:   if (part->vertex_weights) {
137:     PetscMalloc1(mat_localsize,&coarse_vertex_weights);
138:     PetscArraycpy(coarse_vertex_weights,part->vertex_weights,mat_localsize);
139:     MatPartitioningSetVertexWeights(hpart->coarseMatPart,coarse_vertex_weights);
140:   }
141:   /* Copy use_edge_weights flag from part to coarse part */
142:   MatPartitioningGetUseEdgeWeights(part,&use_edge_weights);
143:   MatPartitioningSetUseEdgeWeights(hpart->coarseMatPart,use_edge_weights);

145:   MatPartitioningSetPartitionWeights(hpart->coarseMatPart, part_weights);
146:   MatPartitioningApply(hpart->coarseMatPart,&hpart->coarseparts);

148:   /* Wrap the original vertex weights into an index set so that we can extract the corresponding
149:    * vertex weights for each big subdomain using ISCreateSubIS().
150:    * */
151:   if (part->vertex_weights) {
152:     ISCreateGeneral(comm,mat_localsize,part->vertex_weights,PETSC_COPY_VALUES,&vweights);
153:   }

155:   PetscCalloc1(mat_localsize, &fineparts_indices_tmp);
156:   for (i=0; i<hpart->ncoarseparts; i+=size) {
157:     /* Determine where we want to send big subdomains */
158:     MatPartitioningHierarchical_DetermineDestination(part,hpart->coarseparts,i,i+size,&destination);
159:     /* Assemble a submatrix and its vertex weights for partitioning subdomains  */
160:     MatPartitioningHierarchical_AssembleSubdomain(adj,part->vertex_weights? vweights:NULL,destination,part->vertex_weights? &svweights:NULL,&sadj,&mapping);
161:     /* We have to create a new array to hold vertex weights since coarse partitioner needs to own the vertex-weights array */
162:     if (part->vertex_weights) {
163:       ISGetLocalSize(svweights,&nsvwegihts);
164:       PetscMalloc1(nsvwegihts,&fp_vweights);
165:       ISGetIndices(svweights,&svweights_indices);
166:       PetscArraycpy(fp_vweights,svweights_indices,nsvwegihts);
167:       ISRestoreIndices(svweights,&svweights_indices);
168:       ISDestroy(&svweights);
169:     }

171:     ISDestroy(&destination);
172:     PetscObjectGetComm((PetscObject)sadj,&scomm);

174:     /*
175:      * If the number of big subdomains is smaller than the number of processor cores, the higher ranks do not
176:      * need to do partitioning
177:      * */
178:     if ((i+rank)<hpart->ncoarseparts) {
179:       MatPartitioningDestroy(&hpart->fineMatPart);
180:       /* create a fine partitioner */
181:       MatPartitioningCreate(scomm,&hpart->fineMatPart);
182:       PetscObjectSetOptionsPrefix((PetscObject)hpart->fineMatPart,prefix);
183:       PetscObjectAppendOptionsPrefix((PetscObject)hpart->fineMatPart,"hierarch_fine_");
184:       /* if do not set partitioning type, use parmetis by default */
185:       if (!hpart->fineparttype) {
186: #if defined(PETSC_HAVE_PARMETIS)
187:         MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPARMETIS);
188:         PetscStrallocpy(MATPARTITIONINGPARMETIS,&hpart->fineparttype);
189: #elif defined(PETSC_HAVE_PTSCOTCH)
190:         MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPTSCOTCH);
191:         PetscStrallocpy(MATPARTITIONINGPTSCOTCH,&hpart->fineparttype);
192: #elif defined(PETSC_HAVE_CHACO)
193:         MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGCHACO);
194:         PetscStrallocpy(MATPARTITIONINGCHACO,&hpart->fineparttype);
195: #elif defined(PETSC_HAVE_PARTY)
196:         MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPARTY);
197:         PetscStrallocpy(PETSC_HAVE_PARTY,&hpart->fineparttype);
198: #else
199:         SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
200: #endif
201:       } else {
202:         MatPartitioningSetType(hpart->fineMatPart,hpart->fineparttype);
203:       }
204:       MatPartitioningSetUseEdgeWeights(hpart->fineMatPart,use_edge_weights);
205:       MatPartitioningSetAdjacency(hpart->fineMatPart,sadj);
206:       MatPartitioningSetNParts(hpart->fineMatPart, offsets[rank+1+i]-offsets[rank+i]);
207:       if (part->vertex_weights) {
208:         MatPartitioningSetVertexWeights(hpart->fineMatPart,fp_vweights);
209:       }
210:       MatPartitioningApply(hpart->fineMatPart,&fineparts_temp);
211:     } else {
212:       ISCreateGeneral(scomm,0,NULL,PETSC_OWN_POINTER,&fineparts_temp);
213:     }

215:     MatDestroy(&sadj);

217:     /* Send partition back to the original owners */
218:     MatPartitioningHierarchical_ReassembleFineparts(adj,fineparts_temp,mapping,&hpart->fineparts);
219:     ISGetIndices(hpart->fineparts,&fineparts_indices);
220:     for (j=0;j<mat_localsize;j++)
221:       if (fineparts_indices[j] >=0) fineparts_indices_tmp[j] = fineparts_indices[j];

223:     ISRestoreIndices(hpart->fineparts,&fineparts_indices);
224:     ISDestroy(&hpart->fineparts);
225:     ISDestroy(&fineparts_temp);
226:     ISLocalToGlobalMappingDestroy(&mapping);
227:   }

229:   if (part->vertex_weights) {
230:     ISDestroy(&vweights);
231:   }

233:   ISCreateGeneral(comm,mat_localsize,fineparts_indices_tmp,PETSC_OWN_POINTER,&hpart->fineparts);
234:   ISGetIndices(hpart->fineparts,&fineparts_indices);
235:   ISGetIndices(hpart->coarseparts,&coarseparts_indices);
236:   PetscMalloc1(bs*adj->rmap->n,&parts_indices);
237:   /* Modify the local indices to the global indices by combing the coarse partition and the fine partitions */
238:   for (i=0; i<adj->rmap->n; i++) {
239:     for (j=0; j<bs; j++) {
240:       parts_indices[bs*i+j] = fineparts_indices[i]+offsets[coarseparts_indices[i]];
241:     }
242:   }
243:   ISRestoreIndices(hpart->fineparts,&fineparts_indices);
244:   ISRestoreIndices(hpart->coarseparts,&coarseparts_indices);
245:   PetscFree(offsets);
246:   ISCreateGeneral(comm,bs*adj->rmap->n,parts_indices,PETSC_OWN_POINTER,partitioning);
247:   MatDestroy(&adj);
248:   return 0;
249: }

251: PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat adj, IS fineparts, ISLocalToGlobalMapping mapping, IS *sfineparts)
252: {
253:   PetscInt            *local_indices, *global_indices,*sfineparts_indices,localsize,i;
254:   const PetscInt      *ranges,*fineparts_indices;
255:   PetscMPIInt         rank,*owners;
256:   MPI_Comm            comm;
257:   PetscLayout         rmap;
258:   PetscSFNode        *remote;
259:   PetscSF             sf;

262:   PetscObjectGetComm((PetscObject)adj,&comm);
263:   MPI_Comm_rank(comm,&rank);
264:   MatGetLayouts(adj,&rmap,NULL);
265:   ISGetLocalSize(fineparts,&localsize);
266:   PetscMalloc2(localsize,&global_indices,localsize,&local_indices);
267:   for (i=0; i<localsize; i++) {
268:     local_indices[i] = i;
269:   }
270:   /* map local indices back to global so that we can permulate data globally */
271:   ISLocalToGlobalMappingApply(mapping,localsize,local_indices,global_indices);
272:   PetscCalloc1(localsize,&owners);
273:   /* find owners for global indices */
274:   for (i=0; i<localsize; i++) {
275:     PetscLayoutFindOwner(rmap,global_indices[i],&owners[i]);
276:   }
277:   PetscLayoutGetRanges(rmap,&ranges);
278:   PetscMalloc1(ranges[rank+1]-ranges[rank],&sfineparts_indices);

280:   for (i=0; i<ranges[rank+1]-ranges[rank]; i++) {
281:     sfineparts_indices[i] = -1;
282:   }

284:   ISGetIndices(fineparts,&fineparts_indices);
285:   PetscSFCreate(comm,&sf);
286:   PetscMalloc1(localsize,&remote);
287:   for (i=0; i<localsize; i++) {
288:     remote[i].rank  = owners[i];
289:     remote[i].index = global_indices[i]-ranges[owners[i]];
290:   }
291:   PetscSFSetType(sf,PETSCSFBASIC);
292:   /* not sure how to add prefix to sf */
293:   PetscSFSetFromOptions(sf);
294:   PetscSFSetGraph(sf,localsize,localsize,NULL,PETSC_OWN_POINTER,remote,PETSC_OWN_POINTER);
295:   PetscSFReduceBegin(sf,MPIU_INT,fineparts_indices,sfineparts_indices,MPI_REPLACE);
296:   PetscSFReduceEnd(sf,MPIU_INT,fineparts_indices,sfineparts_indices,MPI_REPLACE);
297:   PetscSFDestroy(&sf);
298:   ISRestoreIndices(fineparts,&fineparts_indices);
299:   ISCreateGeneral(comm,ranges[rank+1]-ranges[rank],sfineparts_indices,PETSC_OWN_POINTER,sfineparts);
300:   PetscFree2(global_indices,local_indices);
301:   PetscFree(owners);
302:   return 0;
303: }

305: PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat adj,IS vweights, IS destination,IS *svweights,Mat *sadj,ISLocalToGlobalMapping *mapping)
306: {
307:   IS              irows,icols;
308:   PetscInt        irows_ln;
309:   PetscMPIInt     rank;
310:   const PetscInt *irows_indices;
311:   MPI_Comm        comm;

313:   PetscObjectGetComm((PetscObject)adj,&comm);
314:   MPI_Comm_rank(comm,&rank);
315:   /* figure out where data comes from  */
316:   ISBuildTwoSided(destination,NULL,&irows);
317:   ISDuplicate(irows,&icols);
318:   ISGetLocalSize(irows,&irows_ln);
319:   ISGetIndices(irows,&irows_indices);
320:   ISLocalToGlobalMappingCreate(comm,1,irows_ln,irows_indices,PETSC_COPY_VALUES,mapping);
321:   ISRestoreIndices(irows,&irows_indices);
322:   MatCreateSubMatrices(adj,1,&irows,&icols,MAT_INITIAL_MATRIX,&sadj);
323:   if (vweights && svweights) {
324:     ISCreateSubIS(vweights,irows,svweights);
325:   }
326:   ISDestroy(&irows);
327:   ISDestroy(&icols);
328:   return 0;
329: }

331: PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning part, IS partitioning, PetscInt pstart, PetscInt pend, IS *destination)
332: {
333:   MPI_Comm            comm;
334:   PetscMPIInt         rank,size,target;
335:   PetscInt            plocalsize,*dest_indices,i;
336:   const PetscInt     *part_indices;

338:   PetscObjectGetComm((PetscObject)part,&comm);
339:   MPI_Comm_rank(comm,&rank);
340:   MPI_Comm_size(comm,&size);
343:   ISGetLocalSize(partitioning,&plocalsize);
344:   PetscMalloc1(plocalsize,&dest_indices);
345:   ISGetIndices(partitioning,&part_indices);
346:   for (i=0; i<plocalsize; i++) {
347:     /* compute target */
348:     target = part_indices[i]-pstart;
349:     /* mark out of range entity as -1 */
350:     if (part_indices[i]<pstart || part_indices[i]>=pend) target = -1;
351:     dest_indices[i] = target;
352:   }
353:   ISCreateGeneral(comm,plocalsize,dest_indices,PETSC_OWN_POINTER,destination);
354:   return 0;
355: }

357: PetscErrorCode MatPartitioningView_Hierarchical(MatPartitioning part,PetscViewer viewer)
358: {
359:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
360:   PetscMPIInt              rank;
361:   PetscBool                iascii;
362:   PetscViewer              sviewer;

364:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
365:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
366:   if (iascii) {
367:     PetscViewerASCIIPrintf(viewer," Number of coarse parts: %" PetscInt_FMT "\n",hpart->ncoarseparts);
368:     PetscViewerASCIIPrintf(viewer," Coarse partitioner: %s\n",hpart->coarseparttype);
369:     if (hpart->coarseMatPart) {
370:       PetscViewerASCIIPushTab(viewer);
371:       MatPartitioningView(hpart->coarseMatPart,viewer);
372:       PetscViewerASCIIPopTab(viewer);
373:     }
374:     PetscViewerASCIIPrintf(viewer," Number of fine parts: %" PetscInt_FMT "\n",hpart->nfineparts);
375:     PetscViewerASCIIPrintf(viewer," Fine partitioner: %s\n",hpart->fineparttype);
376:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
377:     if (rank == 0 && hpart->fineMatPart) {
378:       PetscViewerASCIIPushTab(viewer);
379:       MatPartitioningView(hpart->fineMatPart,sviewer);
380:       PetscViewerASCIIPopTab(viewer);
381:     }
382:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
383:   }
384:   return 0;
385: }

387: PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning part,IS *fineparts)
388: {
389:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

391:   *fineparts = hpart->fineparts;
392:   PetscObjectReference((PetscObject)hpart->fineparts);
393:   return 0;
394: }

396: PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning part,IS *coarseparts)
397: {
398:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

400:   *coarseparts = hpart->coarseparts;
401:   PetscObjectReference((PetscObject)hpart->coarseparts);
402:   return 0;
403: }

405: PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning part, PetscInt ncoarseparts)
406: {
407:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

409:   hpart->ncoarseparts = ncoarseparts;
410:   return 0;
411: }

413: PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning part, PetscInt nfineparts)
414: {
415:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

417:   hpart->nfineparts = nfineparts;
418:   return 0;
419: }

421: PetscErrorCode MatPartitioningSetFromOptions_Hierarchical(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
422: {
423:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
424:   char           value[1024];
425:   PetscBool      flag = PETSC_FALSE;

427:   PetscOptionsHead(PetscOptionsObject,"Set hierarchical partitioning options");
428:   PetscOptionsString("-mat_partitioning_hierarchical_coarseparttype","coarse part type",NULL,NULL,value,sizeof(value),&flag);
429:   if (flag) {
430:    PetscStrallocpy(value,&hpart->coarseparttype);
431:   }
432:   PetscOptionsString("-mat_partitioning_hierarchical_fineparttype","fine part type",NULL,NULL,value,sizeof(value),&flag);
433:   if (flag) {
434:     PetscStrallocpy(value,&hpart->fineparttype);
435:   }
436:   PetscOptionsInt("-mat_partitioning_hierarchical_ncoarseparts","number of coarse parts",NULL,hpart->ncoarseparts,&hpart->ncoarseparts,&flag);
437:   PetscOptionsInt("-mat_partitioning_hierarchical_nfineparts","number of fine parts",NULL,hpart->nfineparts,&hpart->nfineparts,&flag);
438:   PetscOptionsTail();
439:   return 0;
440: }

442: PetscErrorCode MatPartitioningDestroy_Hierarchical(MatPartitioning part)
443: {
444:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

446:   if (hpart->coarseparttype) PetscFree(hpart->coarseparttype);
447:   if (hpart->fineparttype) PetscFree(hpart->fineparttype);
448:   ISDestroy(&hpart->fineparts);
449:   ISDestroy(&hpart->coarseparts);
450:   MatPartitioningDestroy(&hpart->coarseMatPart);
451:   MatPartitioningDestroy(&hpart->fineMatPart);
452:   MatPartitioningDestroy(&hpart->improver);
453:   PetscFree(hpart);
454:   return 0;
455: }

457: /*
458:    Improves the quality  of a partition
459: */
460: static PetscErrorCode MatPartitioningImprove_Hierarchical(MatPartitioning part, IS *partitioning)
461: {
462:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
463:   Mat                           mat = part->adj, adj;
464:   PetscBool                    flg;
465:   const char                   *prefix;
466: #if defined(PETSC_HAVE_PARMETIS)
467:   PetscInt                     *vertex_weights;
468: #endif

470:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
471:   if (flg) {
472:     adj = mat;
473:     PetscObjectReference((PetscObject)adj);
474:   }else {
475:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
476:        resulting partition results need to be stretched to match the original matrix */
477:    MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
478:   }

480:   /* If there exists a mat partitioner, we should delete it */
481:   MatPartitioningDestroy(&hpart->improver);
482:   MatPartitioningCreate(PetscObjectComm((PetscObject)part),&hpart->improver);
483:   PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
484:   PetscObjectSetOptionsPrefix((PetscObject)hpart->improver,prefix);
485:   PetscObjectAppendOptionsPrefix((PetscObject)hpart->improver,"hierarch_improver_");
486:   /* Only parmetis supports to refine a partition */
487: #if defined(PETSC_HAVE_PARMETIS)
488:   MatPartitioningSetType(hpart->improver,MATPARTITIONINGPARMETIS);
489:   MatPartitioningSetAdjacency(hpart->improver,adj);
490:   MatPartitioningSetNParts(hpart->improver, part->n);
491:   /* copy over vertex weights */
492:   if (part->vertex_weights) {
493:     PetscMalloc1(adj->rmap->n,&vertex_weights);
494:     PetscArraycpy(vertex_weights,part->vertex_weights,adj->rmap->n);
495:     MatPartitioningSetVertexWeights(hpart->improver,vertex_weights);
496:   }
497:   MatPartitioningImprove(hpart->improver,partitioning);
498:   MatDestroy(&adj);
499:   return 0;
500: #else
501:   SETERRQ(PetscObjectComm((PetscObject)adj),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis");
502: #endif
503: }

505: /*MC
506:    MATPARTITIONINGHIERARCH - Creates a partitioning context via hierarchical partitioning strategy.
507:    The graph is partitioned into a number of subgraphs, and each subgraph is further split into a few smaller
508:    subgraphs. The idea can be applied in a recursive manner. It is useful when you want to partition the graph
509:    into a large number of subgraphs (often more than 10K) since partitions obtained with existing partitioners
510:    such as ParMETIS and PTScotch are far from ideal. The hierarchical partitioning also tries to avoid off-node
511:    communication as much as possible for multi-core processor. Another user case for the hierarchical partitioning
512:    is to improve PCGASM convergence by generating multi-rank connected subdomain.

514:    Collective

516:    Input Parameter:
517: .  part - the partitioning context

519:    Options Database Keys:
520: +     -mat_partitioning_hierarchical_coarseparttype - partitioner type at the first level and parmetis is used by default
521: .     -mat_partitioning_hierarchical_fineparttype - partitioner type at the second level and parmetis is used by default
522: .     -mat_partitioning_hierarchical_ncoarseparts - number of subgraphs is required at the first level, which is often the number of compute nodes
523: -     -mat_partitioning_hierarchical_nfineparts - number of smaller subgraphs for each subgraph, which is often the number of cores per compute node

525:    Level: beginner

527:    References:
528: +  * - Fande Kong, Xiao-Chuan Cai, A highly scalable multilevel Schwarz method with boundary geometry preserving coarse spaces for 3D elasticity
529:       problems on domains with complex geometry,   SIAM Journal on Scientific Computing 38 (2), C73-C95, 2016
530: -  * - Fande Kong, Roy H. Stogner, Derek Gaston, John W. Peterson, Cody J. Permann, Andrew E. Slaughter, and Richard C. Martineau,
531:       A general-purpose hierarchical mesh partitioning method with node balancing strategies for large-scale numerical simulations,
532:       arXiv preprint arXiv:1809.02666CoRR, 2018.

534: .seealso: MatPartitioningSetType(), MatPartitioningType

536: M*/

538: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning part)
539: {
540:   MatPartitioning_Hierarchical *hpart;

542:   PetscNewLog(part,&hpart);
543:   part->data = (void*)hpart;

545:   hpart->fineparttype       = NULL; /* fine level (second) partitioner */
546:   hpart->coarseparttype     = NULL; /* coarse level (first) partitioner */
547:   hpart->nfineparts         = 1;    /* we do not further partition coarse partition any more by default */
548:   hpart->ncoarseparts       = 0;    /* number of coarse parts (first level) */
549:   hpart->coarseparts        = NULL;
550:   hpart->fineparts          = NULL;
551:   hpart->coarseMatPart      = NULL;
552:   hpart->fineMatPart        = NULL;

554:   part->ops->apply          = MatPartitioningApply_Hierarchical;
555:   part->ops->view           = MatPartitioningView_Hierarchical;
556:   part->ops->destroy        = MatPartitioningDestroy_Hierarchical;
557:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Hierarchical;
558:   part->ops->improve        = MatPartitioningImprove_Hierarchical;
559:   return 0;
560: }