Actual source code: hierarchical.c
petsc-3.12.2 2019-11-22
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 = 0;
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: PetscErrorCode ierr;
56: PetscObjectGetComm((PetscObject)part,&comm);
57: MPI_Comm_size(comm,&size);
58: MPI_Comm_rank(comm,&rank);
59: PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
60: if (flg) {
61: adj = mat;
62: PetscObjectReference((PetscObject)adj);
63: }else {
64: /* bs indicates if the converted matrix is "reduced" from the original and hence the
65: resulting partition results need to be stretched to match the original matrix */
66: MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
67: if (adj->rmap->n > 0) bs = mat->rmap->n/adj->rmap->n;
68: }
69: /* local size of mat */
70: mat_localsize = adj->rmap->n;
71: /* check parameters */
72: /* how many small subdomains we want from a given 'big' suddomain */
73: if(!hpart->nfineparts) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," must set number of small subdomains for each big subdomain \n");
74: 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");
76: /* Partitioning the domain into one single subdomain is a trivial case, and we should just return */
77: if (part->n==1) {
78: PetscCalloc1(bs*adj->rmap->n,&parts_indices);
79: ISCreateGeneral(comm,bs*adj->rmap->n,parts_indices,PETSC_OWN_POINTER,partitioning);
80: hpart->ncoarseparts = 1;
81: hpart->nfineparts = 1;
82: PetscStrallocpy("NONE",&hpart->coarseparttype);
83: PetscStrallocpy("NONE",&hpart->fineparttype);
84: MatDestroy(&adj);
85: return(0);
86: }
88: if(part->n){
89: hpart->ncoarseparts = part->n/hpart->nfineparts;
91: if (part->n%hpart->nfineparts != 0) hpart->ncoarseparts++;
92: }else{
93: part->n = hpart->ncoarseparts*hpart->nfineparts;
94: }
96: PetscMalloc1(hpart->ncoarseparts+1, &offsets);
97: PetscMalloc1(hpart->ncoarseparts, &part_weights);
99: offsets[0] = 0;
100: if (part->n%hpart->nfineparts != 0) offsets[1] = part->n%hpart->nfineparts;
101: else offsets[1] = hpart->nfineparts;
103: part_weights[0] = ((PetscReal)offsets[1])/part->n;
105: for (i=2; i<=hpart->ncoarseparts; i++) {
106: offsets[i] = hpart->nfineparts;
107: part_weights[i-1] = ((PetscReal)offsets[i])/part->n;
108: }
110: offsets[0] = 0;
111: for (i=1;i<=hpart->ncoarseparts; i++)
112: offsets[i] += offsets[i-1];
114: /* If these exists a mat partitioner, we should delete it */
115: MatPartitioningDestroy(&hpart->coarseMatPart);
116: MatPartitioningCreate(comm,&hpart->coarseMatPart);
117: PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
118: PetscObjectSetOptionsPrefix((PetscObject)hpart->coarseMatPart,prefix);
119: PetscObjectAppendOptionsPrefix((PetscObject)hpart->coarseMatPart,"hierarch_coarse_");
120: /* if did not set partitioning type yet, use parmetis by default */
121: if (!hpart->coarseparttype){
122: #if defined(PETSC_HAVE_PARMETIS)
123: MatPartitioningSetType(hpart->coarseMatPart,MATPARTITIONINGPARMETIS);
124: PetscStrallocpy(MATPARTITIONINGPARMETIS,&hpart->coarseparttype);
125: #elif defined(PETSC_HAVE_PTSCOTCH)
126: MatPartitioningSetType(hpart->coarseMatPart,MATPARTITIONINGPTSCOTCH);
127: PetscStrallocpy(MATPARTITIONINGPTSCOTCH,&hpart->coarseparttype);
128: #else
129: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
130: #endif
131: } else {
132: MatPartitioningSetType(hpart->coarseMatPart,hpart->coarseparttype);
133: }
134: MatPartitioningSetAdjacency(hpart->coarseMatPart,adj);
135: MatPartitioningSetNParts(hpart->coarseMatPart, hpart->ncoarseparts);
136: /* copy over vertex weights */
137: if(part->vertex_weights){
138: PetscMalloc1(mat_localsize,&coarse_vertex_weights);
139: PetscArraycpy(coarse_vertex_weights,part->vertex_weights,mat_localsize);
140: MatPartitioningSetVertexWeights(hpart->coarseMatPart,coarse_vertex_weights);
141: }
143: MatPartitioningSetPartitionWeights(hpart->coarseMatPart, part_weights);
144: MatPartitioningApply(hpart->coarseMatPart,&hpart->coarseparts);
146: /* Wrap the original vertex weights into an index set so that we can extract the corresponding
147: * vertex weights for each big subdomain using ISCreateSubIS().
148: * */
149: if (part->vertex_weights) {
150: ISCreateGeneral(comm,mat_localsize,part->vertex_weights,PETSC_COPY_VALUES,&vweights);
151: }
153: PetscCalloc1(mat_localsize, &fineparts_indices_tmp);
154: for(i=0; i<hpart->ncoarseparts; i+=size){
155: /* Determine where we want to send big subdomains */
156: MatPartitioningHierarchical_DetermineDestination(part,hpart->coarseparts,i,i+size,&destination);
157: /* Assemble a submatrix and its vertex weights for partitioning subdomains */
158: MatPartitioningHierarchical_AssembleSubdomain(adj,part->vertex_weights? vweights:NULL,destination,part->vertex_weights? &svweights:NULL,&sadj,&mapping);
159: /* We have to create a new array to hold vertex weights since coarse partitioner needs to own the vertex-weights array */
160: if (part->vertex_weights) {
161: ISGetLocalSize(svweights,&nsvwegihts);
162: PetscMalloc1(nsvwegihts,&fp_vweights);
163: ISGetIndices(svweights,&svweights_indices);
164: PetscArraycpy(fp_vweights,svweights_indices,nsvwegihts);
165: ISRestoreIndices(svweights,&svweights_indices);
166: ISDestroy(&svweights);
167: }
169: ISDestroy(&destination);
170: PetscObjectGetComm((PetscObject)sadj,&scomm);
172: /*
173: * If the number of big subdomains is smaller than the number of processor cores, the higher ranks do not
174: * need to do partitioning
175: * */
176: if((i+rank)<hpart->ncoarseparts) {
177: MatPartitioningDestroy(&hpart->fineMatPart);
178: /* create a fine partitioner */
179: MatPartitioningCreate(scomm,&hpart->fineMatPart);
180: PetscObjectSetOptionsPrefix((PetscObject)hpart->fineMatPart,prefix);
181: PetscObjectAppendOptionsPrefix((PetscObject)hpart->fineMatPart,"hierarch_fine_");
182: /* if do not set partitioning type, use parmetis by default */
183: if(!hpart->fineparttype){
184: #if defined(PETSC_HAVE_PARMETIS)
185: MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPARMETIS);
186: PetscStrallocpy(MATPARTITIONINGPARMETIS,&hpart->fineparttype);
187: #elif defined(PETSC_HAVE_PTSCOTCH)
188: MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPTSCOTCH);
189: PetscStrallocpy(MATPARTITIONINGPTSCOTCH,&hpart->fineparttype);
190: #elif defined(PETSC_HAVE_CHACO)
191: MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGCHACO);
192: PetscStrallocpy(MATPARTITIONINGCHACO,&hpart->fineparttype);
193: #elif defined(PETSC_HAVE_PARTY)
194: MatPartitioningSetType(hpart->fineMatPart,MATPARTITIONINGPARTY);
195: PetscStrallocpy(PETSC_HAVE_PARTY,&hpart->fineparttype);
196: #else
197: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
198: #endif
199: } else {
200: MatPartitioningSetType(hpart->fineMatPart,hpart->fineparttype);
201: }
202: MatPartitioningSetAdjacency(hpart->fineMatPart,sadj);
203: MatPartitioningSetNParts(hpart->fineMatPart, offsets[rank+1+i]-offsets[rank+i]);
204: if (part->vertex_weights) {
205: MatPartitioningSetVertexWeights(hpart->fineMatPart,fp_vweights);
206: }
207: MatPartitioningApply(hpart->fineMatPart,&fineparts_temp);
208: } else {
209: ISCreateGeneral(scomm,0,NULL,PETSC_OWN_POINTER,&fineparts_temp);
210: }
212: MatDestroy(&sadj);
214: /* Send partition back to the original owners */
215: MatPartitioningHierarchical_ReassembleFineparts(adj,fineparts_temp,mapping,&hpart->fineparts);
216: ISGetIndices(hpart->fineparts,&fineparts_indices);
217: for (j=0;j<mat_localsize;j++)
218: if (fineparts_indices[j] >=0) fineparts_indices_tmp[j] = fineparts_indices[j];
220: ISRestoreIndices(hpart->fineparts,&fineparts_indices);
221: ISDestroy(&hpart->fineparts);
222: ISDestroy(&fineparts_temp);
223: ISLocalToGlobalMappingDestroy(&mapping);
224: }
226: if (part->vertex_weights) {
227: ISDestroy(&vweights);
228: }
230: ISCreateGeneral(comm,mat_localsize,fineparts_indices_tmp,PETSC_OWN_POINTER,&hpart->fineparts);
231: ISGetIndices(hpart->fineparts,&fineparts_indices);
232: ISGetIndices(hpart->coarseparts,&coarseparts_indices);
233: PetscMalloc1(bs*adj->rmap->n,&parts_indices);
234: /* Modify the local indices to the global indices by combing the coarse partition and the fine partitions */
235: for(i=0; i<adj->rmap->n; i++){
236: for(j=0; j<bs; j++){
237: parts_indices[bs*i+j] = fineparts_indices[i]+offsets[coarseparts_indices[i]];
238: }
239: }
240: ISRestoreIndices(hpart->fineparts,&fineparts_indices);
241: ISRestoreIndices(hpart->coarseparts,&coarseparts_indices);
242: PetscFree(offsets);
243: ISCreateGeneral(comm,bs*adj->rmap->n,parts_indices,PETSC_OWN_POINTER,partitioning);
244: MatDestroy(&adj);
245: return(0);
246: }
249: PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat adj, IS fineparts, ISLocalToGlobalMapping mapping, IS *sfineparts)
250: {
251: PetscInt *local_indices, *global_indices,*owners,*sfineparts_indices,localsize,i;
252: const PetscInt *ranges,*fineparts_indices;
253: PetscMPIInt rank;
254: MPI_Comm comm;
255: PetscLayout rmap;
256: PetscSFNode *remote;
257: PetscSF sf;
258: PetscErrorCode ierr;
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,MPIU_REPLACE);
296: PetscSFReduceEnd(sf,MPIU_INT,fineparts_indices,sfineparts_indices,MPIU_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: }
306: PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat adj,IS vweights, IS destination,IS *svweights,Mat *sadj,ISLocalToGlobalMapping *mapping)
307: {
308: IS irows,icols;
309: PetscInt irows_ln;
310: PetscMPIInt rank;
311: const PetscInt *irows_indices;
312: MPI_Comm comm;
313: PetscErrorCode ierr;
316: PetscObjectGetComm((PetscObject)adj,&comm);
317: MPI_Comm_rank(comm,&rank);
318: /* figure out where data comes from */
319: ISBuildTwoSided(destination,NULL,&irows);
320: ISDuplicate(irows,&icols);
321: ISGetLocalSize(irows,&irows_ln);
322: ISGetIndices(irows,&irows_indices);
323: ISLocalToGlobalMappingCreate(comm,1,irows_ln,irows_indices,PETSC_COPY_VALUES,mapping);
324: ISRestoreIndices(irows,&irows_indices);
325: MatCreateSubMatrices(adj,1,&irows,&icols,MAT_INITIAL_MATRIX,&sadj);
326: if (vweights && svweights) {
327: ISCreateSubIS(vweights,irows,svweights);
328: }
329: ISDestroy(&irows);
330: ISDestroy(&icols);
331: return(0);
332: }
335: PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning part, IS partitioning, PetscInt pstart, PetscInt pend, IS *destination)
336: {
337: MPI_Comm comm;
338: PetscMPIInt rank,size,target;
339: PetscInt plocalsize,*dest_indices,i;
340: const PetscInt *part_indices;
341: PetscErrorCode ierr;
344: PetscObjectGetComm((PetscObject)part,&comm);
345: MPI_Comm_rank(comm,&rank);
346: MPI_Comm_size(comm,&size);
347: 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);
348: if(pstart>pend) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," pstart %D should be smaller than pend %D",pstart,pend);
349: ISGetLocalSize(partitioning,&plocalsize);
350: PetscMalloc1(plocalsize,&dest_indices);
351: ISGetIndices(partitioning,&part_indices);
352: for (i=0; i<plocalsize; i++){
353: /* compute target */
354: target = part_indices[i]-pstart;
355: /* mark out of range entity as -1 */
356: if (part_indices[i]<pstart || part_indices[i]>=pend) target = -1;
357: dest_indices[i] = target;
358: }
359: ISCreateGeneral(comm,plocalsize,dest_indices,PETSC_OWN_POINTER,destination);
360: return(0);
361: }
364: PetscErrorCode MatPartitioningView_Hierarchical(MatPartitioning part,PetscViewer viewer)
365: {
366: MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
367: PetscErrorCode ierr;
368: PetscMPIInt rank;
369: PetscBool iascii;
370: PetscViewer sviewer;
373: MPI_Comm_rank(PetscObjectComm((PetscObject)part),&rank);
374: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
375: if(iascii){
376: PetscViewerASCIIPrintf(viewer," Number of coarse parts: %D\n",hpart->ncoarseparts);
377: PetscViewerASCIIPrintf(viewer," Coarse partitioner: %s\n",hpart->coarseparttype);
378: if (hpart->coarseMatPart) {
379: PetscViewerASCIIPushTab(viewer);
380: MatPartitioningView(hpart->coarseMatPart,viewer);
381: PetscViewerASCIIPopTab(viewer);
382: }
383: PetscViewerASCIIPrintf(viewer," Number of fine parts: %D\n",hpart->nfineparts);
384: PetscViewerASCIIPrintf(viewer," Fine partitioner: %s\n",hpart->fineparttype);
385: PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
386: if (!rank && hpart->fineMatPart) {
387: PetscViewerASCIIPushTab(viewer);
388: MatPartitioningView(hpart->fineMatPart,sviewer);
389: PetscViewerASCIIPopTab(viewer);
390: }
391: PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);
392: }
393: return(0);
394: }
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,1024,&flag);
447: if (flag){
448: PetscStrallocpy(value,&hpart->coarseparttype);
449: }
450: PetscOptionsString("-mat_partitioning_hierarchical_fineparttype","fine part type",NULL,NULL,value,1024,&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: }
461: PetscErrorCode MatPartitioningDestroy_Hierarchical(MatPartitioning part)
462: {
463: MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
464: PetscErrorCode ierr;
467: if(hpart->coarseparttype) {PetscFree(hpart->coarseparttype);}
468: if(hpart->fineparttype) {PetscFree(hpart->fineparttype);}
469: ISDestroy(&hpart->fineparts);
470: ISDestroy(&hpart->coarseparts);
471: MatPartitioningDestroy(&hpart->coarseMatPart);
472: MatPartitioningDestroy(&hpart->fineMatPart);
473: MatPartitioningDestroy(&hpart->improver);
474: PetscFree(hpart);
475: return(0);
476: }
478: /*
479: Improves the quality of a partition
480: */
481: static PetscErrorCode MatPartitioningImprove_Hierarchical(MatPartitioning part, IS *partitioning)
482: {
483: PetscErrorCode ierr;
484: MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical*)part->data;
485: Mat mat = part->adj, adj;
486: PetscBool flg;
487: PetscInt *vertex_weights;
488: const char *prefix;
491: PetscObjectTypeCompare((PetscObject)mat,MATMPIADJ,&flg);
492: if (flg) {
493: adj = mat;
494: PetscObjectReference((PetscObject)adj);
495: }else {
496: /* bs indicates if the converted matrix is "reduced" from the original and hence the
497: resulting partition results need to be stretched to match the original matrix */
498: MatConvert(mat,MATMPIADJ,MAT_INITIAL_MATRIX,&adj);
499: }
501: /* If there exists a mat partitioner, we should delete it */
502: MatPartitioningDestroy(&hpart->improver);
503: MatPartitioningCreate(PetscObjectComm((PetscObject)part),&hpart->improver);
504: PetscObjectGetOptionsPrefix((PetscObject)part,&prefix);
505: PetscObjectSetOptionsPrefix((PetscObject)hpart->improver,prefix);
506: PetscObjectAppendOptionsPrefix((PetscObject)hpart->improver,"hierarch_improver_");
507: /* Only parmetis supports to refine a partition */
508: #if defined(PETSC_HAVE_PARMETIS)
509: MatPartitioningSetType(hpart->improver,MATPARTITIONINGPARMETIS);
510: #else
511: SETERRQ(PetscObjectComm((PetscObject)adj),PETSC_ERR_SUP,"Requires PETSc be installed with ParMetis\n");
512: #endif
514: MatPartitioningSetAdjacency(hpart->improver,adj);
515: MatPartitioningSetNParts(hpart->improver, part->n);
516: /* copy over vertex weights */
517: if(part->vertex_weights){
518: PetscMalloc1(adj->rmap->n,&vertex_weights);
519: PetscArraycpy(vertex_weights,part->vertex_weights,adj->rmap->n);
520: MatPartitioningSetVertexWeights(hpart->improver,vertex_weights);
521: }
522: MatPartitioningImprove(hpart->improver,partitioning);
523: MatDestroy(&adj);
524: return(0);
525: }
528: /*MC
529: MATPARTITIONINGHIERARCH - Creates a partitioning context via hierarchical partitioning strategy.
530: The graph is partitioned into a number of subgraphs, and each subgraph is further split into a few smaller
531: subgraphs. The idea can be applied in a recursive manner. It is useful when you want to partition the graph
532: into a large number of subgraphs (often more than 10K) since partitions obtained with existing partitioners
533: such as ParMETIS and PTScotch are far from ideal. The hierarchical partitioning also tries to avoid off-node
534: communication as much as possible for multi-core processor. Another user case for the hierarchical partitioning
535: is to improve PCGASM convergence by generating multi-rank connected subdomain.
537: Collective
539: Input Parameter:
540: . part - the partitioning context
542: Options Database Keys:
543: + -mat_partitioning_hierarchical_coarseparttype - partitioner type at the first level and parmetis is used by default
544: . -mat_partitioning_hierarchical_fineparttype - partitioner type at the second level and parmetis is used by default
545: . -mat_partitioning_hierarchical_ncoarseparts - number of subgraphs is required at the first level, which is often the number of compute nodes
546: - -mat_partitioning_hierarchical_nfineparts - number of smaller subgraphs for each subgraph, which is often the number of cores per compute node
548: Level: beginner
550: References:
551: 1. Fande Kong, Xiao-Chuan Cai, A highly scalable multilevel Schwarz method with boundary geometry preserving coarse spaces for 3D elasticity
552: problems on domains with complex geometry, SIAM Journal on Scientific Computing 38 (2), C73-C95, 2016
553: 2. Fande Kong, Roy H. Stogner, Derek Gaston, John W. Peterson, Cody J. Permann, Andrew E. Slaughter, and Richard C. Martineau,
554: A general-purpose hierarchical mesh partitioning method with node balancing strategies for large-scale numerical simulations,
555: arXiv preprint arXiv:1809.02666CoRR, 2018.
557: .seealso: MatPartitioningSetType(), MatPartitioningType
559: M*/
561: PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning part)
562: {
563: PetscErrorCode ierr;
564: MatPartitioning_Hierarchical *hpart;
567: PetscNewLog(part,&hpart);
568: part->data = (void*)hpart;
570: hpart->fineparttype = 0; /* fine level (second) partitioner */
571: hpart->coarseparttype = 0; /* coarse level (first) partitioner */
572: hpart->nfineparts = 1; /* we do not further partition coarse partition any more by default */
573: hpart->ncoarseparts = 0; /* number of coarse parts (first level) */
574: hpart->coarseparts = 0;
575: hpart->fineparts = 0;
576: hpart->coarseMatPart = 0;
577: hpart->fineMatPart = 0;
579: part->ops->apply = MatPartitioningApply_Hierarchical;
580: part->ops->view = MatPartitioningView_Hierarchical;
581: part->ops->destroy = MatPartitioningDestroy_Hierarchical;
582: part->ops->setfromoptions = MatPartitioningSetFromOptions_Hierarchical;
583: part->ops->improve = MatPartitioningImprove_Hierarchical;
584: return(0);
585: }