xref: /petsc/src/mat/graphops/partition/impls/hierarchical/hierarchical.c (revision 8be712e46db5d855f641c6bd97b4543e0efe65bd)
1*8be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
2*8be712e4SBarry Smith #include <petscsf.h>
3*8be712e4SBarry Smith #include <petsc/private/matimpl.h>
4*8be712e4SBarry Smith 
5*8be712e4SBarry Smith /*
6*8be712e4SBarry Smith   It is a hierarchical partitioning. The partitioner has two goals:
7*8be712e4SBarry Smith   (1) Most of current partitioners fail at a large scale. The hierarchical partitioning
8*8be712e4SBarry Smith   strategy is trying to produce large number of subdomains when number of processor cores is large.
9*8be712e4SBarry Smith   (2) PCGASM needs one 'big' subdomain across multi-cores. The partitioner provides two
10*8be712e4SBarry Smith   consistent partitions, coarse parts and fine parts. A coarse part is a 'big' subdomain consisting
11*8be712e4SBarry Smith   of several small subdomains.
12*8be712e4SBarry Smith */
13*8be712e4SBarry Smith 
14*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning, IS, PetscInt, PetscInt, IS *);
15*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat, IS, IS, IS *, Mat *, ISLocalToGlobalMapping *);
16*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat, IS, ISLocalToGlobalMapping, IS *);
17*8be712e4SBarry Smith 
18*8be712e4SBarry Smith typedef struct {
19*8be712e4SBarry Smith   char           *fineparttype;   /* partitioner on fine level */
20*8be712e4SBarry Smith   char           *coarseparttype; /* partitioner on coarse level */
21*8be712e4SBarry Smith   PetscInt        nfineparts;     /* number of fine parts on each coarse subdomain */
22*8be712e4SBarry Smith   PetscInt        ncoarseparts;   /* number of coarse parts */
23*8be712e4SBarry Smith   IS              coarseparts;    /* partitioning on coarse level */
24*8be712e4SBarry Smith   IS              fineparts;      /* partitioning on fine level */
25*8be712e4SBarry Smith   MatPartitioning coarseMatPart;  /* MatPartititioning on coarse level (first level) */
26*8be712e4SBarry Smith   MatPartitioning fineMatPart;    /* MatPartitioning on fine level (second level) */
27*8be712e4SBarry Smith   MatPartitioning improver;       /* Improve the quality of a partition */
28*8be712e4SBarry Smith } MatPartitioning_Hierarchical;
29*8be712e4SBarry Smith 
30*8be712e4SBarry Smith /*
31*8be712e4SBarry Smith    Uses a hierarchical partitioning strategy to partition the matrix in parallel.
32*8be712e4SBarry Smith    Use this interface to make the partitioner consistent with others
33*8be712e4SBarry Smith */
34*8be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Hierarchical(MatPartitioning part, IS *partitioning)
35*8be712e4SBarry Smith {
36*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
37*8be712e4SBarry Smith   const PetscInt               *fineparts_indices, *coarseparts_indices;
38*8be712e4SBarry Smith   PetscInt                     *fineparts_indices_tmp;
39*8be712e4SBarry Smith   PetscInt                     *parts_indices, i, j, mat_localsize, *offsets;
40*8be712e4SBarry Smith   Mat                           mat = part->adj, adj, sadj;
41*8be712e4SBarry Smith   PetscReal                    *part_weights;
42*8be712e4SBarry Smith   PetscBool                     flg;
43*8be712e4SBarry Smith   PetscInt                      bs                    = 1;
44*8be712e4SBarry Smith   PetscInt                     *coarse_vertex_weights = NULL;
45*8be712e4SBarry Smith   PetscMPIInt                   size, rank;
46*8be712e4SBarry Smith   MPI_Comm                      comm, scomm;
47*8be712e4SBarry Smith   IS                            destination, fineparts_temp, vweights, svweights;
48*8be712e4SBarry Smith   PetscInt                      nsvwegihts, *fp_vweights;
49*8be712e4SBarry Smith   const PetscInt               *svweights_indices;
50*8be712e4SBarry Smith   ISLocalToGlobalMapping        mapping;
51*8be712e4SBarry Smith   const char                   *prefix;
52*8be712e4SBarry Smith   PetscBool                     use_edge_weights;
53*8be712e4SBarry Smith 
54*8be712e4SBarry Smith   PetscFunctionBegin;
55*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
56*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
57*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
58*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
59*8be712e4SBarry Smith   if (flg) {
60*8be712e4SBarry Smith     adj = mat;
61*8be712e4SBarry Smith     PetscCall(PetscObjectReference((PetscObject)adj));
62*8be712e4SBarry Smith   } else {
63*8be712e4SBarry Smith     /* bs indicates if the converted matrix is "reduced" from the original and hence the
64*8be712e4SBarry Smith        resulting partition results need to be stretched to match the original matrix */
65*8be712e4SBarry Smith     PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
66*8be712e4SBarry Smith     if (adj->rmap->n > 0) bs = mat->rmap->n / adj->rmap->n;
67*8be712e4SBarry Smith   }
68*8be712e4SBarry Smith   /* local size of mat */
69*8be712e4SBarry Smith   mat_localsize = adj->rmap->n;
70*8be712e4SBarry Smith   /* check parameters */
71*8be712e4SBarry Smith   /* how many small subdomains we want from a given 'big' suddomain */
72*8be712e4SBarry Smith   PetscCheck(hpart->nfineparts, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, " must set number of small subdomains for each big subdomain ");
73*8be712e4SBarry Smith   PetscCheck(hpart->ncoarseparts || part->n, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, " did not either set number of coarse parts or total number of parts ");
74*8be712e4SBarry Smith 
75*8be712e4SBarry Smith   /* Partitioning the domain into one single subdomain is a trivial case, and we should just return  */
76*8be712e4SBarry Smith   if (part->n == 1) {
77*8be712e4SBarry Smith     PetscCall(PetscCalloc1(bs * adj->rmap->n, &parts_indices));
78*8be712e4SBarry Smith     PetscCall(ISCreateGeneral(comm, bs * adj->rmap->n, parts_indices, PETSC_OWN_POINTER, partitioning));
79*8be712e4SBarry Smith     hpart->ncoarseparts = 1;
80*8be712e4SBarry Smith     hpart->nfineparts   = 1;
81*8be712e4SBarry Smith     PetscCall(PetscStrallocpy("NONE", &hpart->coarseparttype));
82*8be712e4SBarry Smith     PetscCall(PetscStrallocpy("NONE", &hpart->fineparttype));
83*8be712e4SBarry Smith     PetscCall(MatDestroy(&adj));
84*8be712e4SBarry Smith     PetscFunctionReturn(PETSC_SUCCESS);
85*8be712e4SBarry Smith   }
86*8be712e4SBarry Smith 
87*8be712e4SBarry Smith   if (part->n) {
88*8be712e4SBarry Smith     hpart->ncoarseparts = part->n / hpart->nfineparts;
89*8be712e4SBarry Smith 
90*8be712e4SBarry Smith     if (part->n % hpart->nfineparts != 0) hpart->ncoarseparts++;
91*8be712e4SBarry Smith   } else {
92*8be712e4SBarry Smith     part->n = hpart->ncoarseparts * hpart->nfineparts;
93*8be712e4SBarry Smith   }
94*8be712e4SBarry Smith 
95*8be712e4SBarry Smith   PetscCall(PetscMalloc1(hpart->ncoarseparts + 1, &offsets));
96*8be712e4SBarry Smith   PetscCall(PetscMalloc1(hpart->ncoarseparts, &part_weights));
97*8be712e4SBarry Smith 
98*8be712e4SBarry Smith   offsets[0] = 0;
99*8be712e4SBarry Smith   if (part->n % hpart->nfineparts != 0) offsets[1] = part->n % hpart->nfineparts;
100*8be712e4SBarry Smith   else offsets[1] = hpart->nfineparts;
101*8be712e4SBarry Smith 
102*8be712e4SBarry Smith   part_weights[0] = ((PetscReal)offsets[1]) / part->n;
103*8be712e4SBarry Smith 
104*8be712e4SBarry Smith   for (i = 2; i <= hpart->ncoarseparts; i++) {
105*8be712e4SBarry Smith     offsets[i]          = hpart->nfineparts;
106*8be712e4SBarry Smith     part_weights[i - 1] = ((PetscReal)offsets[i]) / part->n;
107*8be712e4SBarry Smith   }
108*8be712e4SBarry Smith 
109*8be712e4SBarry Smith   offsets[0] = 0;
110*8be712e4SBarry Smith   for (i = 1; i <= hpart->ncoarseparts; i++) offsets[i] += offsets[i - 1];
111*8be712e4SBarry Smith 
112*8be712e4SBarry Smith   /* If these exists a mat partitioner, we should delete it */
113*8be712e4SBarry Smith   PetscCall(MatPartitioningDestroy(&hpart->coarseMatPart));
114*8be712e4SBarry Smith   PetscCall(MatPartitioningCreate(comm, &hpart->coarseMatPart));
115*8be712e4SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)part, &prefix));
116*8be712e4SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->coarseMatPart, prefix));
117*8be712e4SBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->coarseMatPart, "hierarch_coarse_"));
118*8be712e4SBarry Smith   /* if did not set partitioning type yet, use parmetis by default */
119*8be712e4SBarry Smith   if (!hpart->coarseparttype) {
120*8be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
121*8be712e4SBarry Smith     PetscCall(MatPartitioningSetType(hpart->coarseMatPart, MATPARTITIONINGPARMETIS));
122*8be712e4SBarry Smith     PetscCall(PetscStrallocpy(MATPARTITIONINGPARMETIS, &hpart->coarseparttype));
123*8be712e4SBarry Smith #elif defined(PETSC_HAVE_PTSCOTCH)
124*8be712e4SBarry Smith     PetscCall(MatPartitioningSetType(hpart->coarseMatPart, MATPARTITIONINGPTSCOTCH));
125*8be712e4SBarry Smith     PetscCall(PetscStrallocpy(MATPARTITIONINGPTSCOTCH, &hpart->coarseparttype));
126*8be712e4SBarry Smith #else
127*8be712e4SBarry Smith     SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
128*8be712e4SBarry Smith #endif
129*8be712e4SBarry Smith   } else {
130*8be712e4SBarry Smith     PetscCall(MatPartitioningSetType(hpart->coarseMatPart, hpart->coarseparttype));
131*8be712e4SBarry Smith   }
132*8be712e4SBarry Smith   PetscCall(MatPartitioningSetAdjacency(hpart->coarseMatPart, adj));
133*8be712e4SBarry Smith   PetscCall(MatPartitioningSetNParts(hpart->coarseMatPart, hpart->ncoarseparts));
134*8be712e4SBarry Smith   /* copy over vertex weights */
135*8be712e4SBarry Smith   if (part->vertex_weights) {
136*8be712e4SBarry Smith     PetscCall(PetscMalloc1(mat_localsize, &coarse_vertex_weights));
137*8be712e4SBarry Smith     PetscCall(PetscArraycpy(coarse_vertex_weights, part->vertex_weights, mat_localsize));
138*8be712e4SBarry Smith     PetscCall(MatPartitioningSetVertexWeights(hpart->coarseMatPart, coarse_vertex_weights));
139*8be712e4SBarry Smith   }
140*8be712e4SBarry Smith   /* Copy use_edge_weights flag from part to coarse part */
141*8be712e4SBarry Smith   PetscCall(MatPartitioningGetUseEdgeWeights(part, &use_edge_weights));
142*8be712e4SBarry Smith   PetscCall(MatPartitioningSetUseEdgeWeights(hpart->coarseMatPart, use_edge_weights));
143*8be712e4SBarry Smith 
144*8be712e4SBarry Smith   PetscCall(MatPartitioningSetPartitionWeights(hpart->coarseMatPart, part_weights));
145*8be712e4SBarry Smith   PetscCall(MatPartitioningApply(hpart->coarseMatPart, &hpart->coarseparts));
146*8be712e4SBarry Smith 
147*8be712e4SBarry Smith   /* Wrap the original vertex weights into an index set so that we can extract the corresponding
148*8be712e4SBarry Smith    * vertex weights for each big subdomain using ISCreateSubIS().
149*8be712e4SBarry Smith    * */
150*8be712e4SBarry Smith   if (part->vertex_weights) PetscCall(ISCreateGeneral(comm, mat_localsize, part->vertex_weights, PETSC_COPY_VALUES, &vweights));
151*8be712e4SBarry Smith 
152*8be712e4SBarry Smith   PetscCall(PetscCalloc1(mat_localsize, &fineparts_indices_tmp));
153*8be712e4SBarry Smith   for (i = 0; i < hpart->ncoarseparts; i += size) {
154*8be712e4SBarry Smith     /* Determine where we want to send big subdomains */
155*8be712e4SBarry Smith     PetscCall(MatPartitioningHierarchical_DetermineDestination(part, hpart->coarseparts, i, i + size, &destination));
156*8be712e4SBarry Smith     /* Assemble a submatrix and its vertex weights for partitioning subdomains  */
157*8be712e4SBarry Smith     PetscCall(MatPartitioningHierarchical_AssembleSubdomain(adj, part->vertex_weights ? vweights : NULL, destination, part->vertex_weights ? &svweights : NULL, &sadj, &mapping));
158*8be712e4SBarry Smith     /* We have to create a new array to hold vertex weights since coarse partitioner needs to own the vertex-weights array */
159*8be712e4SBarry Smith     if (part->vertex_weights) {
160*8be712e4SBarry Smith       PetscCall(ISGetLocalSize(svweights, &nsvwegihts));
161*8be712e4SBarry Smith       PetscCall(PetscMalloc1(nsvwegihts, &fp_vweights));
162*8be712e4SBarry Smith       PetscCall(ISGetIndices(svweights, &svweights_indices));
163*8be712e4SBarry Smith       PetscCall(PetscArraycpy(fp_vweights, svweights_indices, nsvwegihts));
164*8be712e4SBarry Smith       PetscCall(ISRestoreIndices(svweights, &svweights_indices));
165*8be712e4SBarry Smith       PetscCall(ISDestroy(&svweights));
166*8be712e4SBarry Smith     }
167*8be712e4SBarry Smith 
168*8be712e4SBarry Smith     PetscCall(ISDestroy(&destination));
169*8be712e4SBarry Smith     PetscCall(PetscObjectGetComm((PetscObject)sadj, &scomm));
170*8be712e4SBarry Smith 
171*8be712e4SBarry Smith     /*
172*8be712e4SBarry Smith      * If the number of big subdomains is smaller than the number of processor cores, the higher ranks do not
173*8be712e4SBarry Smith      * need to do partitioning
174*8be712e4SBarry Smith      * */
175*8be712e4SBarry Smith     if ((i + rank) < hpart->ncoarseparts) {
176*8be712e4SBarry Smith       PetscCall(MatPartitioningDestroy(&hpart->fineMatPart));
177*8be712e4SBarry Smith       /* create a fine partitioner */
178*8be712e4SBarry Smith       PetscCall(MatPartitioningCreate(scomm, &hpart->fineMatPart));
179*8be712e4SBarry Smith       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->fineMatPart, prefix));
180*8be712e4SBarry Smith       PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->fineMatPart, "hierarch_fine_"));
181*8be712e4SBarry Smith       /* if do not set partitioning type, use parmetis by default */
182*8be712e4SBarry Smith       if (!hpart->fineparttype) {
183*8be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
184*8be712e4SBarry Smith         PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPARMETIS));
185*8be712e4SBarry Smith         PetscCall(PetscStrallocpy(MATPARTITIONINGPARMETIS, &hpart->fineparttype));
186*8be712e4SBarry Smith #elif defined(PETSC_HAVE_PTSCOTCH)
187*8be712e4SBarry Smith         PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPTSCOTCH));
188*8be712e4SBarry Smith         PetscCall(PetscStrallocpy(MATPARTITIONINGPTSCOTCH, &hpart->fineparttype));
189*8be712e4SBarry Smith #elif defined(PETSC_HAVE_CHACO)
190*8be712e4SBarry Smith         PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGCHACO));
191*8be712e4SBarry Smith         PetscCall(PetscStrallocpy(MATPARTITIONINGCHACO, &hpart->fineparttype));
192*8be712e4SBarry Smith #elif defined(PETSC_HAVE_PARTY)
193*8be712e4SBarry Smith         PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPARTY));
194*8be712e4SBarry Smith         PetscCall(PetscStrallocpy(PETSC_HAVE_PARTY, &hpart->fineparttype));
195*8be712e4SBarry Smith #else
196*8be712e4SBarry Smith         SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
197*8be712e4SBarry Smith #endif
198*8be712e4SBarry Smith       } else {
199*8be712e4SBarry Smith         PetscCall(MatPartitioningSetType(hpart->fineMatPart, hpart->fineparttype));
200*8be712e4SBarry Smith       }
201*8be712e4SBarry Smith       PetscCall(MatPartitioningSetUseEdgeWeights(hpart->fineMatPart, use_edge_weights));
202*8be712e4SBarry Smith       PetscCall(MatPartitioningSetAdjacency(hpart->fineMatPart, sadj));
203*8be712e4SBarry Smith       PetscCall(MatPartitioningSetNParts(hpart->fineMatPart, offsets[rank + 1 + i] - offsets[rank + i]));
204*8be712e4SBarry Smith       if (part->vertex_weights) PetscCall(MatPartitioningSetVertexWeights(hpart->fineMatPart, fp_vweights));
205*8be712e4SBarry Smith       PetscCall(MatPartitioningApply(hpart->fineMatPart, &fineparts_temp));
206*8be712e4SBarry Smith     } else {
207*8be712e4SBarry Smith       PetscCall(ISCreateGeneral(scomm, 0, NULL, PETSC_OWN_POINTER, &fineparts_temp));
208*8be712e4SBarry Smith     }
209*8be712e4SBarry Smith 
210*8be712e4SBarry Smith     PetscCall(MatDestroy(&sadj));
211*8be712e4SBarry Smith 
212*8be712e4SBarry Smith     /* Send partition back to the original owners */
213*8be712e4SBarry Smith     PetscCall(MatPartitioningHierarchical_ReassembleFineparts(adj, fineparts_temp, mapping, &hpart->fineparts));
214*8be712e4SBarry Smith     PetscCall(ISGetIndices(hpart->fineparts, &fineparts_indices));
215*8be712e4SBarry Smith     for (j = 0; j < mat_localsize; j++)
216*8be712e4SBarry Smith       if (fineparts_indices[j] >= 0) fineparts_indices_tmp[j] = fineparts_indices[j];
217*8be712e4SBarry Smith 
218*8be712e4SBarry Smith     PetscCall(ISRestoreIndices(hpart->fineparts, &fineparts_indices));
219*8be712e4SBarry Smith     PetscCall(ISDestroy(&hpart->fineparts));
220*8be712e4SBarry Smith     PetscCall(ISDestroy(&fineparts_temp));
221*8be712e4SBarry Smith     PetscCall(ISLocalToGlobalMappingDestroy(&mapping));
222*8be712e4SBarry Smith   }
223*8be712e4SBarry Smith 
224*8be712e4SBarry Smith   if (part->vertex_weights) PetscCall(ISDestroy(&vweights));
225*8be712e4SBarry Smith 
226*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(comm, mat_localsize, fineparts_indices_tmp, PETSC_OWN_POINTER, &hpart->fineparts));
227*8be712e4SBarry Smith   PetscCall(ISGetIndices(hpart->fineparts, &fineparts_indices));
228*8be712e4SBarry Smith   PetscCall(ISGetIndices(hpart->coarseparts, &coarseparts_indices));
229*8be712e4SBarry Smith   PetscCall(PetscMalloc1(bs * adj->rmap->n, &parts_indices));
230*8be712e4SBarry Smith   /* Modify the local indices to the global indices by combing the coarse partition and the fine partitions */
231*8be712e4SBarry Smith   for (i = 0; i < adj->rmap->n; i++) {
232*8be712e4SBarry Smith     for (j = 0; j < bs; j++) parts_indices[bs * i + j] = fineparts_indices[i] + offsets[coarseparts_indices[i]];
233*8be712e4SBarry Smith   }
234*8be712e4SBarry Smith   PetscCall(ISRestoreIndices(hpart->fineparts, &fineparts_indices));
235*8be712e4SBarry Smith   PetscCall(ISRestoreIndices(hpart->coarseparts, &coarseparts_indices));
236*8be712e4SBarry Smith   PetscCall(PetscFree(offsets));
237*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(comm, bs * adj->rmap->n, parts_indices, PETSC_OWN_POINTER, partitioning));
238*8be712e4SBarry Smith   PetscCall(MatDestroy(&adj));
239*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
240*8be712e4SBarry Smith }
241*8be712e4SBarry Smith 
242*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat adj, IS fineparts, ISLocalToGlobalMapping mapping, IS *sfineparts)
243*8be712e4SBarry Smith {
244*8be712e4SBarry Smith   PetscInt       *local_indices, *global_indices, *sfineparts_indices, localsize, i;
245*8be712e4SBarry Smith   const PetscInt *ranges, *fineparts_indices;
246*8be712e4SBarry Smith   PetscMPIInt     rank, *owners;
247*8be712e4SBarry Smith   MPI_Comm        comm;
248*8be712e4SBarry Smith   PetscLayout     rmap;
249*8be712e4SBarry Smith   PetscSFNode    *remote;
250*8be712e4SBarry Smith   PetscSF         sf;
251*8be712e4SBarry Smith 
252*8be712e4SBarry Smith   PetscFunctionBegin;
253*8be712e4SBarry Smith   PetscAssertPointer(sfineparts, 4);
254*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)adj, &comm));
255*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
256*8be712e4SBarry Smith   PetscCall(MatGetLayouts(adj, &rmap, NULL));
257*8be712e4SBarry Smith   PetscCall(ISGetLocalSize(fineparts, &localsize));
258*8be712e4SBarry Smith   PetscCall(PetscMalloc2(localsize, &global_indices, localsize, &local_indices));
259*8be712e4SBarry Smith   for (i = 0; i < localsize; i++) local_indices[i] = i;
260*8be712e4SBarry Smith   /* map local indices back to global so that we can permulate data globally */
261*8be712e4SBarry Smith   PetscCall(ISLocalToGlobalMappingApply(mapping, localsize, local_indices, global_indices));
262*8be712e4SBarry Smith   PetscCall(PetscCalloc1(localsize, &owners));
263*8be712e4SBarry Smith   /* find owners for global indices */
264*8be712e4SBarry Smith   for (i = 0; i < localsize; i++) PetscCall(PetscLayoutFindOwner(rmap, global_indices[i], &owners[i]));
265*8be712e4SBarry Smith   PetscCall(PetscLayoutGetRanges(rmap, &ranges));
266*8be712e4SBarry Smith   PetscCall(PetscMalloc1(ranges[rank + 1] - ranges[rank], &sfineparts_indices));
267*8be712e4SBarry Smith 
268*8be712e4SBarry Smith   for (i = 0; i < ranges[rank + 1] - ranges[rank]; i++) sfineparts_indices[i] = -1;
269*8be712e4SBarry Smith 
270*8be712e4SBarry Smith   PetscCall(ISGetIndices(fineparts, &fineparts_indices));
271*8be712e4SBarry Smith   PetscCall(PetscSFCreate(comm, &sf));
272*8be712e4SBarry Smith   PetscCall(PetscMalloc1(localsize, &remote));
273*8be712e4SBarry Smith   for (i = 0; i < localsize; i++) {
274*8be712e4SBarry Smith     remote[i].rank  = owners[i];
275*8be712e4SBarry Smith     remote[i].index = global_indices[i] - ranges[owners[i]];
276*8be712e4SBarry Smith   }
277*8be712e4SBarry Smith   PetscCall(PetscSFSetType(sf, PETSCSFBASIC));
278*8be712e4SBarry Smith   /* not sure how to add prefix to sf */
279*8be712e4SBarry Smith   PetscCall(PetscSFSetFromOptions(sf));
280*8be712e4SBarry Smith   PetscCall(PetscSFSetGraph(sf, localsize, localsize, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
281*8be712e4SBarry Smith   PetscCall(PetscSFReduceBegin(sf, MPIU_INT, fineparts_indices, sfineparts_indices, MPI_REPLACE));
282*8be712e4SBarry Smith   PetscCall(PetscSFReduceEnd(sf, MPIU_INT, fineparts_indices, sfineparts_indices, MPI_REPLACE));
283*8be712e4SBarry Smith   PetscCall(PetscSFDestroy(&sf));
284*8be712e4SBarry Smith   PetscCall(ISRestoreIndices(fineparts, &fineparts_indices));
285*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(comm, ranges[rank + 1] - ranges[rank], sfineparts_indices, PETSC_OWN_POINTER, sfineparts));
286*8be712e4SBarry Smith   PetscCall(PetscFree2(global_indices, local_indices));
287*8be712e4SBarry Smith   PetscCall(PetscFree(owners));
288*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
289*8be712e4SBarry Smith }
290*8be712e4SBarry Smith 
291*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat adj, IS vweights, IS destination, IS *svweights, Mat *sadj, ISLocalToGlobalMapping *mapping)
292*8be712e4SBarry Smith {
293*8be712e4SBarry Smith   IS              irows, icols;
294*8be712e4SBarry Smith   PetscInt        irows_ln;
295*8be712e4SBarry Smith   PetscMPIInt     rank;
296*8be712e4SBarry Smith   const PetscInt *irows_indices;
297*8be712e4SBarry Smith   MPI_Comm        comm;
298*8be712e4SBarry Smith 
299*8be712e4SBarry Smith   PetscFunctionBegin;
300*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)adj, &comm));
301*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
302*8be712e4SBarry Smith   /* figure out where data comes from  */
303*8be712e4SBarry Smith   PetscCall(ISBuildTwoSided(destination, NULL, &irows));
304*8be712e4SBarry Smith   PetscCall(ISDuplicate(irows, &icols));
305*8be712e4SBarry Smith   PetscCall(ISGetLocalSize(irows, &irows_ln));
306*8be712e4SBarry Smith   PetscCall(ISGetIndices(irows, &irows_indices));
307*8be712e4SBarry Smith   PetscCall(ISLocalToGlobalMappingCreate(comm, 1, irows_ln, irows_indices, PETSC_COPY_VALUES, mapping));
308*8be712e4SBarry Smith   PetscCall(ISRestoreIndices(irows, &irows_indices));
309*8be712e4SBarry Smith   PetscCall(MatCreateSubMatrices(adj, 1, &irows, &icols, MAT_INITIAL_MATRIX, &sadj));
310*8be712e4SBarry Smith   if (vweights && svweights) PetscCall(ISCreateSubIS(vweights, irows, svweights));
311*8be712e4SBarry Smith   PetscCall(ISDestroy(&irows));
312*8be712e4SBarry Smith   PetscCall(ISDestroy(&icols));
313*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
314*8be712e4SBarry Smith }
315*8be712e4SBarry Smith 
316*8be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning part, IS partitioning, PetscInt pstart, PetscInt pend, IS *destination)
317*8be712e4SBarry Smith {
318*8be712e4SBarry Smith   MPI_Comm        comm;
319*8be712e4SBarry Smith   PetscMPIInt     rank, size, target;
320*8be712e4SBarry Smith   PetscInt        plocalsize, *dest_indices, i;
321*8be712e4SBarry Smith   const PetscInt *part_indices;
322*8be712e4SBarry Smith 
323*8be712e4SBarry Smith   PetscFunctionBegin;
324*8be712e4SBarry Smith   PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
325*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(comm, &rank));
326*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_size(comm, &size));
327*8be712e4SBarry Smith   PetscCheck((pend - pstart) <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "range [%" PetscInt_FMT ", %" PetscInt_FMT "] should be smaller than or equal to size %d", pstart, pend, size);
328*8be712e4SBarry Smith   PetscCheck(pstart <= pend, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, " pstart %" PetscInt_FMT " should be smaller than pend %" PetscInt_FMT, pstart, pend);
329*8be712e4SBarry Smith   PetscCall(ISGetLocalSize(partitioning, &plocalsize));
330*8be712e4SBarry Smith   PetscCall(PetscMalloc1(plocalsize, &dest_indices));
331*8be712e4SBarry Smith   PetscCall(ISGetIndices(partitioning, &part_indices));
332*8be712e4SBarry Smith   for (i = 0; i < plocalsize; i++) {
333*8be712e4SBarry Smith     /* compute target */
334*8be712e4SBarry Smith     target = part_indices[i] - pstart;
335*8be712e4SBarry Smith     /* mark out of range entity as -1 */
336*8be712e4SBarry Smith     if (part_indices[i] < pstart || part_indices[i] >= pend) target = -1;
337*8be712e4SBarry Smith     dest_indices[i] = target;
338*8be712e4SBarry Smith   }
339*8be712e4SBarry Smith   PetscCall(ISCreateGeneral(comm, plocalsize, dest_indices, PETSC_OWN_POINTER, destination));
340*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
341*8be712e4SBarry Smith }
342*8be712e4SBarry Smith 
343*8be712e4SBarry Smith static PetscErrorCode MatPartitioningView_Hierarchical(MatPartitioning part, PetscViewer viewer)
344*8be712e4SBarry Smith {
345*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
346*8be712e4SBarry Smith   PetscMPIInt                   rank;
347*8be712e4SBarry Smith   PetscBool                     iascii;
348*8be712e4SBarry Smith   PetscViewer                   sviewer;
349*8be712e4SBarry Smith 
350*8be712e4SBarry Smith   PetscFunctionBegin;
351*8be712e4SBarry Smith   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank));
352*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
353*8be712e4SBarry Smith   if (iascii) {
354*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, " Number of coarse parts: %" PetscInt_FMT "\n", hpart->ncoarseparts));
355*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse partitioner: %s\n", hpart->coarseparttype));
356*8be712e4SBarry Smith     if (hpart->coarseMatPart) {
357*8be712e4SBarry Smith       PetscCall(PetscViewerASCIIPushTab(viewer));
358*8be712e4SBarry Smith       PetscCall(MatPartitioningView(hpart->coarseMatPart, viewer));
359*8be712e4SBarry Smith       PetscCall(PetscViewerASCIIPopTab(viewer));
360*8be712e4SBarry Smith     }
361*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, " Number of fine parts: %" PetscInt_FMT "\n", hpart->nfineparts));
362*8be712e4SBarry Smith     PetscCall(PetscViewerASCIIPrintf(viewer, " Fine partitioner: %s\n", hpart->fineparttype));
363*8be712e4SBarry Smith     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
364*8be712e4SBarry Smith     if (rank == 0 && hpart->fineMatPart) {
365*8be712e4SBarry Smith       PetscCall(PetscViewerASCIIPushTab(viewer));
366*8be712e4SBarry Smith       PetscCall(MatPartitioningView(hpart->fineMatPart, sviewer));
367*8be712e4SBarry Smith       PetscCall(PetscViewerASCIIPopTab(viewer));
368*8be712e4SBarry Smith     }
369*8be712e4SBarry Smith     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
370*8be712e4SBarry Smith   }
371*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
372*8be712e4SBarry Smith }
373*8be712e4SBarry Smith 
374*8be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning part, IS *fineparts)
375*8be712e4SBarry Smith {
376*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
377*8be712e4SBarry Smith 
378*8be712e4SBarry Smith   PetscFunctionBegin;
379*8be712e4SBarry Smith   *fineparts = hpart->fineparts;
380*8be712e4SBarry Smith   PetscCall(PetscObjectReference((PetscObject)hpart->fineparts));
381*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
382*8be712e4SBarry Smith }
383*8be712e4SBarry Smith 
384*8be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning part, IS *coarseparts)
385*8be712e4SBarry Smith {
386*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
387*8be712e4SBarry Smith 
388*8be712e4SBarry Smith   PetscFunctionBegin;
389*8be712e4SBarry Smith   *coarseparts = hpart->coarseparts;
390*8be712e4SBarry Smith   PetscCall(PetscObjectReference((PetscObject)hpart->coarseparts));
391*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
392*8be712e4SBarry Smith }
393*8be712e4SBarry Smith 
394*8be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning part, PetscInt ncoarseparts)
395*8be712e4SBarry Smith {
396*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
397*8be712e4SBarry Smith 
398*8be712e4SBarry Smith   PetscFunctionBegin;
399*8be712e4SBarry Smith   hpart->ncoarseparts = ncoarseparts;
400*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
401*8be712e4SBarry Smith }
402*8be712e4SBarry Smith 
403*8be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning part, PetscInt nfineparts)
404*8be712e4SBarry Smith {
405*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
406*8be712e4SBarry Smith 
407*8be712e4SBarry Smith   PetscFunctionBegin;
408*8be712e4SBarry Smith   hpart->nfineparts = nfineparts;
409*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
410*8be712e4SBarry Smith }
411*8be712e4SBarry Smith 
412*8be712e4SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_Hierarchical(MatPartitioning part, PetscOptionItems *PetscOptionsObject)
413*8be712e4SBarry Smith {
414*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
415*8be712e4SBarry Smith   char                          value[1024];
416*8be712e4SBarry Smith   PetscBool                     flag = PETSC_FALSE;
417*8be712e4SBarry Smith 
418*8be712e4SBarry Smith   PetscFunctionBegin;
419*8be712e4SBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Set hierarchical partitioning options");
420*8be712e4SBarry Smith   PetscCall(PetscOptionsString("-mat_partitioning_hierarchical_coarseparttype", "coarse part type", NULL, NULL, value, sizeof(value), &flag));
421*8be712e4SBarry Smith   if (flag) {
422*8be712e4SBarry Smith     PetscCall(PetscFree(hpart->coarseparttype));
423*8be712e4SBarry Smith     PetscCall(PetscStrallocpy(value, &hpart->coarseparttype));
424*8be712e4SBarry Smith   }
425*8be712e4SBarry Smith   PetscCall(PetscOptionsString("-mat_partitioning_hierarchical_fineparttype", "fine part type", NULL, NULL, value, sizeof(value), &flag));
426*8be712e4SBarry Smith   if (flag) {
427*8be712e4SBarry Smith     PetscCall(PetscFree(hpart->fineparttype));
428*8be712e4SBarry Smith     PetscCall(PetscStrallocpy(value, &hpart->fineparttype));
429*8be712e4SBarry Smith   }
430*8be712e4SBarry Smith   PetscCall(PetscOptionsInt("-mat_partitioning_hierarchical_ncoarseparts", "number of coarse parts", NULL, hpart->ncoarseparts, &hpart->ncoarseparts, &flag));
431*8be712e4SBarry Smith   PetscCall(PetscOptionsInt("-mat_partitioning_hierarchical_nfineparts", "number of fine parts", NULL, hpart->nfineparts, &hpart->nfineparts, &flag));
432*8be712e4SBarry Smith   PetscOptionsHeadEnd();
433*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
434*8be712e4SBarry Smith }
435*8be712e4SBarry Smith 
436*8be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_Hierarchical(MatPartitioning part)
437*8be712e4SBarry Smith {
438*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
439*8be712e4SBarry Smith 
440*8be712e4SBarry Smith   PetscFunctionBegin;
441*8be712e4SBarry Smith   if (hpart->coarseparttype) PetscCall(PetscFree(hpart->coarseparttype));
442*8be712e4SBarry Smith   if (hpart->fineparttype) PetscCall(PetscFree(hpart->fineparttype));
443*8be712e4SBarry Smith   PetscCall(ISDestroy(&hpart->fineparts));
444*8be712e4SBarry Smith   PetscCall(ISDestroy(&hpart->coarseparts));
445*8be712e4SBarry Smith   PetscCall(MatPartitioningDestroy(&hpart->coarseMatPart));
446*8be712e4SBarry Smith   PetscCall(MatPartitioningDestroy(&hpart->fineMatPart));
447*8be712e4SBarry Smith   PetscCall(MatPartitioningDestroy(&hpart->improver));
448*8be712e4SBarry Smith   PetscCall(PetscFree(hpart));
449*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
450*8be712e4SBarry Smith }
451*8be712e4SBarry Smith 
452*8be712e4SBarry Smith /*
453*8be712e4SBarry Smith    Improves the quality  of a partition
454*8be712e4SBarry Smith */
455*8be712e4SBarry Smith static PetscErrorCode MatPartitioningImprove_Hierarchical(MatPartitioning part, IS *partitioning)
456*8be712e4SBarry Smith {
457*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
458*8be712e4SBarry Smith   Mat                           mat   = part->adj, adj;
459*8be712e4SBarry Smith   PetscBool                     flg;
460*8be712e4SBarry Smith   const char                   *prefix;
461*8be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
462*8be712e4SBarry Smith   PetscInt *vertex_weights;
463*8be712e4SBarry Smith #endif
464*8be712e4SBarry Smith 
465*8be712e4SBarry Smith   PetscFunctionBegin;
466*8be712e4SBarry Smith   PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
467*8be712e4SBarry Smith   if (flg) {
468*8be712e4SBarry Smith     adj = mat;
469*8be712e4SBarry Smith     PetscCall(PetscObjectReference((PetscObject)adj));
470*8be712e4SBarry Smith   } else {
471*8be712e4SBarry Smith     /* bs indicates if the converted matrix is "reduced" from the original and hence the
472*8be712e4SBarry Smith        resulting partition results need to be stretched to match the original matrix */
473*8be712e4SBarry Smith     PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
474*8be712e4SBarry Smith   }
475*8be712e4SBarry Smith 
476*8be712e4SBarry Smith   /* If there exists a mat partitioner, we should delete it */
477*8be712e4SBarry Smith   PetscCall(MatPartitioningDestroy(&hpart->improver));
478*8be712e4SBarry Smith   PetscCall(MatPartitioningCreate(PetscObjectComm((PetscObject)part), &hpart->improver));
479*8be712e4SBarry Smith   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)part, &prefix));
480*8be712e4SBarry Smith   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->improver, prefix));
481*8be712e4SBarry Smith   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->improver, "hierarch_improver_"));
482*8be712e4SBarry Smith   /* Only parmetis supports to refine a partition */
483*8be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
484*8be712e4SBarry Smith   PetscCall(MatPartitioningSetType(hpart->improver, MATPARTITIONINGPARMETIS));
485*8be712e4SBarry Smith   PetscCall(MatPartitioningSetAdjacency(hpart->improver, adj));
486*8be712e4SBarry Smith   PetscCall(MatPartitioningSetNParts(hpart->improver, part->n));
487*8be712e4SBarry Smith   /* copy over vertex weights */
488*8be712e4SBarry Smith   if (part->vertex_weights) {
489*8be712e4SBarry Smith     PetscCall(PetscMalloc1(adj->rmap->n, &vertex_weights));
490*8be712e4SBarry Smith     PetscCall(PetscArraycpy(vertex_weights, part->vertex_weights, adj->rmap->n));
491*8be712e4SBarry Smith     PetscCall(MatPartitioningSetVertexWeights(hpart->improver, vertex_weights));
492*8be712e4SBarry Smith   }
493*8be712e4SBarry Smith   PetscCall(MatPartitioningImprove(hpart->improver, partitioning));
494*8be712e4SBarry Smith   PetscCall(MatDestroy(&adj));
495*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
496*8be712e4SBarry Smith #else
497*8be712e4SBarry Smith   SETERRQ(PetscObjectComm((PetscObject)adj), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis");
498*8be712e4SBarry Smith #endif
499*8be712e4SBarry Smith }
500*8be712e4SBarry Smith 
501*8be712e4SBarry Smith /*MC
502*8be712e4SBarry Smith    MATPARTITIONINGHIERARCH - Creates a partitioning context via hierarchical partitioning strategy.
503*8be712e4SBarry Smith    The graph is partitioned into a number of subgraphs, and each subgraph is further split into a few smaller
504*8be712e4SBarry Smith    subgraphs. The idea can be applied in a recursive manner. It is useful when you want to partition the graph
505*8be712e4SBarry Smith    into a large number of subgraphs (often more than 10K) since partitions obtained with existing partitioners
506*8be712e4SBarry Smith    such as ParMETIS and PTScotch are far from ideal. The hierarchical partitioning also tries to avoid off-node
507*8be712e4SBarry Smith    communication as much as possible for multi-core processor. Another user case for the hierarchical partitioning
508*8be712e4SBarry Smith    is to improve `PCGASM` convergence by generating multi-rank connected subdomain.
509*8be712e4SBarry Smith 
510*8be712e4SBarry Smith    Collective
511*8be712e4SBarry Smith 
512*8be712e4SBarry Smith    Input Parameter:
513*8be712e4SBarry Smith .  part - the partitioning context
514*8be712e4SBarry Smith 
515*8be712e4SBarry Smith    Options Database Keys:
516*8be712e4SBarry Smith +     -mat_partitioning_hierarchical_coarseparttype - partitioner type at the first level and parmetis is used by default
517*8be712e4SBarry Smith .     -mat_partitioning_hierarchical_fineparttype - partitioner type at the second level and parmetis is used by default
518*8be712e4SBarry Smith .     -mat_partitioning_hierarchical_ncoarseparts - number of subgraphs is required at the first level, which is often the number of compute nodes
519*8be712e4SBarry Smith -     -mat_partitioning_hierarchical_nfineparts - number of smaller subgraphs for each subgraph, which is often the number of cores per compute node
520*8be712e4SBarry Smith 
521*8be712e4SBarry Smith    Level: beginner
522*8be712e4SBarry Smith 
523*8be712e4SBarry Smith    References:
524*8be712e4SBarry Smith +  * - Fande Kong, Xiao-Chuan Cai, A highly scalable multilevel Schwarz method with boundary geometry preserving coarse spaces for 3D elasticity
525*8be712e4SBarry Smith       problems on domains with complex geometry,   SIAM Journal on Scientific Computing 38 (2), C73-C95, 2016
526*8be712e4SBarry Smith -  * - Fande Kong, Roy H. Stogner, Derek Gaston, John W. Peterson, Cody J. Permann, Andrew E. Slaughter, and Richard C. Martineau,
527*8be712e4SBarry Smith       A general-purpose hierarchical mesh partitioning method with node balancing strategies for large-scale numerical simulations,
528*8be712e4SBarry Smith       arXiv preprint arXiv:1809.02666CoRR, 2018.
529*8be712e4SBarry Smith 
530*8be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MATPARTITIONINGMETIS`, `MATPARTITIONINGPARMETIS`,
531*8be712e4SBarry Smith M*/
532*8be712e4SBarry Smith 
533*8be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning part)
534*8be712e4SBarry Smith {
535*8be712e4SBarry Smith   MatPartitioning_Hierarchical *hpart;
536*8be712e4SBarry Smith 
537*8be712e4SBarry Smith   PetscFunctionBegin;
538*8be712e4SBarry Smith   PetscCall(PetscNew(&hpart));
539*8be712e4SBarry Smith   part->data = (void *)hpart;
540*8be712e4SBarry Smith 
541*8be712e4SBarry Smith   hpart->fineparttype   = NULL; /* fine level (second) partitioner */
542*8be712e4SBarry Smith   hpart->coarseparttype = NULL; /* coarse level (first) partitioner */
543*8be712e4SBarry Smith   hpart->nfineparts     = 1;    /* we do not further partition coarse partition any more by default */
544*8be712e4SBarry Smith   hpart->ncoarseparts   = 0;    /* number of coarse parts (first level) */
545*8be712e4SBarry Smith   hpart->coarseparts    = NULL;
546*8be712e4SBarry Smith   hpart->fineparts      = NULL;
547*8be712e4SBarry Smith   hpart->coarseMatPart  = NULL;
548*8be712e4SBarry Smith   hpart->fineMatPart    = NULL;
549*8be712e4SBarry Smith 
550*8be712e4SBarry Smith   part->ops->apply          = MatPartitioningApply_Hierarchical;
551*8be712e4SBarry Smith   part->ops->view           = MatPartitioningView_Hierarchical;
552*8be712e4SBarry Smith   part->ops->destroy        = MatPartitioningDestroy_Hierarchical;
553*8be712e4SBarry Smith   part->ops->setfromoptions = MatPartitioningSetFromOptions_Hierarchical;
554*8be712e4SBarry Smith   part->ops->improve        = MatPartitioningImprove_Hierarchical;
555*8be712e4SBarry Smith   PetscFunctionReturn(PETSC_SUCCESS);
556*8be712e4SBarry Smith }
557