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;
 54:   PetscErrorCode                ierr;

 57:   PetscObjectGetComm((PetscObject)part,&comm);
 58:   MPI_Comm_size(comm,&size);
 59:   MPI_Comm_rank(comm,&rank);
 60:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
 61:   if (flg) {
 62:     adj = mat;
 63:     PetscObjectReference((PetscObject)adj);
 64:   }else {
 65:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
 66:        resulting partition results need to be stretched to match the original matrix */
 67:    MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
 68:    if (adj->rmap->n > 0) bs = mat->rmap->n/adj->rmap->n;
 69:   }
 70:   /* local size of mat */
 71:   mat_localsize = adj->rmap->n;
 72:   /* check parameters */
 73:   /* how many small subdomains we want from a given 'big' suddomain */
 74:   if (!hpart->nfineparts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," must set number of small subdomains for each big subdomain \n");
 75:   if (!hpart->ncoarseparts && !part->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE," did not either set number of coarse parts or total number of parts \n");

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

 89:   if (part->n) {
 90:     hpart->ncoarseparts = part->n/hpart->nfineparts;

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

 97:   PetscMalloc1(hpart->ncoarseparts+1, &offsets);
 98:   PetscMalloc1(hpart->ncoarseparts, &part_weights);

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

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

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

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

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

147:   MatPartitioningSetPartitionWeights(hpart->coarseMatPart, part_weights);
148:   MatPartitioningApply(hpart->coarseMatPart,&hpart->coarseparts);

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

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

173:     ISDestroy(&destination);
174:     PetscObjectGetComm((PetscObject)sadj,&scomm);

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

217:     MatDestroy(&sadj);

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

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

231:   if (part->vertex_weights) {
232:     ISDestroy(&vweights);
233:   }

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

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

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

284:   for (i=0; i<ranges[rank+1]-ranges[rank]; i++) {
285:     sfineparts_indices[i] = -1;
286:   }

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

309: PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat adj,IS vweights, IS destination,IS *svweights,Mat *sadj,ISLocalToGlobalMapping *mapping)
310: {
311:   IS              irows,icols;
312:   PetscInt        irows_ln;
313:   PetscMPIInt     rank;
314:   const PetscInt *irows_indices;
315:   MPI_Comm        comm;
316:   PetscErrorCode  ierr;

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

337: PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning part, IS partitioning, PetscInt pstart, PetscInt pend, IS *destination)
338: {
339:   MPI_Comm            comm;
340:   PetscMPIInt         rank,size,target;
341:   PetscInt            plocalsize,*dest_indices,i;
342:   const PetscInt     *part_indices;
343:   PetscErrorCode      ierr;

346:   PetscObjectGetComm((PetscObject)part,&comm);
347:   MPI_Comm_rank(comm,&rank);
348:   MPI_Comm_size(comm,&size);
349:   if ((pend-pstart)>size) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"range [%D, %D] should be smaller than or equal to size %D",pstart,pend,size);
350:   if (pstart>pend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," pstart %D should be smaller than pend %D",pstart,pend);
351:   ISGetLocalSize(partitioning,&plocalsize);
352:   PetscMalloc1(plocalsize,&dest_indices);
353:   ISGetIndices(partitioning,&part_indices);
354:   for (i=0; i<plocalsize; i++) {
355:     /* compute target */
356:     target = part_indices[i]-pstart;
357:     /* mark out of range entity as -1 */
358:     if (part_indices[i]<pstart || part_indices[i]>=pend) target = -1;
359:     dest_indices[i] = target;
360:   }
361:   ISCreateGeneral(comm,plocalsize,dest_indices,PETSC_OWN_POINTER,destination);
362:   return(0);
363: }

365: PetscErrorCode MatPartitioningView_Hierarchical(MatPartitioning part,PetscViewer viewer)
366: {
367:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
368:   PetscErrorCode           ierr;
369:   PetscMPIInt              rank;
370:   PetscBool                iascii;
371:   PetscViewer              sviewer;

374:   MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
375:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
376:   if (iascii) {
377:     PetscViewerASCIIPrintf(viewer," Number of coarse parts: %D\n",hpart->ncoarseparts);
378:     PetscViewerASCIIPrintf(viewer," Coarse partitioner: %s\n",hpart->coarseparttype);
379:     if (hpart->coarseMatPart) {
380:       PetscViewerASCIIPushTab(viewer);
381:       MatPartitioningView(hpart->coarseMatPart,viewer);
382:       PetscViewerASCIIPopTab(viewer);
383:     }
384:     PetscViewerASCIIPrintf(viewer," Number of fine parts: %D\n",hpart->nfineparts);
385:     PetscViewerASCIIPrintf(viewer," Fine partitioner: %s\n",hpart->fineparttype);
386:     PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
387:     if (rank == 0 && hpart->fineMatPart) {
388:       PetscViewerASCIIPushTab(viewer);
389:       MatPartitioningView(hpart->fineMatPart,sviewer);
390:       PetscViewerASCIIPopTab(viewer);
391:     }
392:     PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
393:   }
394:   return(0);
395: }

397: PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning part,IS *fineparts)
398: {
399:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
400:   PetscErrorCode                ierr;

403:   *fineparts = hpart->fineparts;
404:   PetscObjectReference((PetscObject)hpart->fineparts);
405:   return(0);
406: }

408: PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning part,IS *coarseparts)
409: {
410:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
411:   PetscErrorCode                ierr;

414:   *coarseparts = hpart->coarseparts;
415:   PetscObjectReference((PetscObject)hpart->coarseparts);
416:   return(0);
417: }

419: PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning part, PetscInt ncoarseparts)
420: {
421:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

424:   hpart->ncoarseparts = ncoarseparts;
425:   return(0);
426: }

428: PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning part, PetscInt nfineparts)
429: {
430:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;

433:   hpart->nfineparts = nfineparts;
434:   return(0);
435: }

437: PetscErrorCode MatPartitioningSetFromOptions_Hierarchical(PetscOptionItems *PetscOptionsObject,MatPartitioning part)
438: {
439:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
441:   char           value[1024];
442:   PetscBool      flag = PETSC_FALSE;

445:   PetscOptionsHead(PetscOptionsObject,"Set hierarchical partitioning options");
446:   PetscOptionsString("-mat_partitioning_hierarchical_coarseparttype","coarse part type",NULL,NULL,value,sizeof(value),&flag);
447:   if (flag) {
448:    PetscStrallocpy(value,&hpart->coarseparttype);
449:   }
450:   PetscOptionsString("-mat_partitioning_hierarchical_fineparttype","fine part type",NULL,NULL,value,sizeof(value),&flag);
451:   if (flag) {
452:     PetscStrallocpy(value,&hpart->fineparttype);
453:   }
454:   PetscOptionsInt("-mat_partitioning_hierarchical_ncoarseparts","number of coarse parts",NULL,hpart->ncoarseparts,&hpart->ncoarseparts,&flag);
455:   PetscOptionsInt("-mat_partitioning_hierarchical_nfineparts","number of fine parts",NULL,hpart->nfineparts,&hpart->nfineparts,&flag);
456:   PetscOptionsTail();
457:   return(0);
458: }

460: PetscErrorCode MatPartitioningDestroy_Hierarchical(MatPartitioning part)
461: {
462:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
463:   PetscErrorCode           ierr;

466:   if (hpart->coarseparttype) {PetscFree(hpart->coarseparttype);}
467:   if (hpart->fineparttype) {PetscFree(hpart->fineparttype);}
468:   ISDestroy(&hpart->fineparts);
469:   ISDestroy(&hpart->coarseparts);
470:   MatPartitioningDestroy(&hpart->coarseMatPart);
471:   MatPartitioningDestroy(&hpart->fineMatPart);
472:   MatPartitioningDestroy(&hpart->improver);
473:   PetscFree(hpart);
474:   return(0);
475: }

477: /*
478:    Improves the quality  of a partition
479: */
480: static PetscErrorCode MatPartitioningImprove_Hierarchical(MatPartitioning part, IS *partitioning)
481: {
482:   PetscErrorCode               ierr;
483:   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
484:   Mat                           mat = part->adj, adj;
485:   PetscBool                    flg;
486:   const char                   *prefix;
487: #if defined(PETSC_HAVE_PARMETIS)
488:   PetscInt                     *vertex_weights;
489: #endif

492:   PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
493:   if (flg) {
494:     adj = mat;
495:     PetscObjectReference((PetscObject)adj);
496:   }else {
497:     /* bs indicates if the converted matrix is "reduced" from the original and hence the
498:        resulting partition results need to be stretched to match the original matrix */
499:    MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
500:   }

502:   /* If there exists a mat partitioner, we should delete it */
503:   MatPartitioningDestroy(&hpart->improver);
504:   MatPartitioningCreate(PetscObjectComm((PetscObject)part),&hpart->improver);
505:   PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
506:   PetscObjectSetOptionsPrefix((PetscObject)hpart->improver,prefix);
507:   PetscObjectAppendOptionsPrefix((PetscObject)hpart->improver,"hierarch_improver_");
508:   /* Only parmetis supports to refine a partition */
509: #if defined(PETSC_HAVE_PARMETIS)
510:   MatPartitioningSetType(hpart->improver,MATPARTITIONINGPARMETIS);
511:   MatPartitioningSetAdjacency(hpart->improver,adj);
512:   MatPartitioningSetNParts(hpart->improver, part->n);
513:   /* copy over vertex weights */
514:   if (part->vertex_weights) {
515:     PetscMalloc1(adj->rmap->n,&vertex_weights);
516:     PetscArraycpy(vertex_weights,part->vertex_weights,adj->rmap->n);
517:     MatPartitioningSetVertexWeights(hpart->improver,vertex_weights);
518:   }
519:   MatPartitioningImprove(hpart->improver,partitioning);
520:   MatDestroy(&adj);
521:   return(0);
522: #else
523:   SETERRQ(PetscObjectComm((PetscObject)adj),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis\n");
524: #endif
525: }

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

536:    Collective

538:    Input Parameter:
539: .  part - the partitioning context

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

547:    Level: beginner

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

556: .seealso: MatPartitioningSetType(), MatPartitioningType

558: M*/

560: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning part)
561: {
562:   PetscErrorCode                ierr;
563:   MatPartitioning_Hierarchical *hpart;

566:   PetscNewLog(part,&hpart);
567:   part->data = (void*)hpart;

569:   hpart->fineparttype       = NULL; /* fine level (second) partitioner */
570:   hpart->coarseparttype     = NULL; /* coarse level (first) partitioner */
571:   hpart->nfineparts         = 1;    /* we do not further partition coarse partition any more by default */
572:   hpart->ncoarseparts       = 0;    /* number of coarse parts (first level) */
573:   hpart->coarseparts        = NULL;
574:   hpart->fineparts          = NULL;
575:   hpart->coarseMatPart      = NULL;
576:   hpart->fineMatPart        = NULL;

578:   part->ops->apply          = MatPartitioningApply_Hierarchical;
579:   part->ops->view           = MatPartitioningView_Hierarchical;
580:   part->ops->destroy        = MatPartitioningDestroy_Hierarchical;
581:   part->ops->setfromoptions = MatPartitioningSetFromOptions_Hierarchical;
582:   part->ops->improve        = MatPartitioningImprove_Hierarchical;
583:   return(0);
584: }