18be712e4SBarry Smith #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/
28be712e4SBarry Smith #include <petscsf.h>
38be712e4SBarry Smith #include <petsc/private/matimpl.h>
48be712e4SBarry Smith
58be712e4SBarry Smith /*
68be712e4SBarry Smith It is a hierarchical partitioning. The partitioner has two goals:
78be712e4SBarry Smith (1) Most of current partitioners fail at a large scale. The hierarchical partitioning
88be712e4SBarry Smith strategy is trying to produce large number of subdomains when number of processor cores is large.
98be712e4SBarry Smith (2) PCGASM needs one 'big' subdomain across multi-cores. The partitioner provides two
108be712e4SBarry Smith consistent partitions, coarse parts and fine parts. A coarse part is a 'big' subdomain consisting
118be712e4SBarry Smith of several small subdomains.
128be712e4SBarry Smith */
138be712e4SBarry Smith
148be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning, IS, PetscInt, PetscInt, IS *);
158be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat, IS, IS, IS *, Mat *, ISLocalToGlobalMapping *);
168be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat, IS, ISLocalToGlobalMapping, IS *);
178be712e4SBarry Smith
188be712e4SBarry Smith typedef struct {
198be712e4SBarry Smith char *fineparttype; /* partitioner on fine level */
208be712e4SBarry Smith char *coarseparttype; /* partitioner on coarse level */
218be712e4SBarry Smith PetscInt nfineparts; /* number of fine parts on each coarse subdomain */
228be712e4SBarry Smith PetscInt ncoarseparts; /* number of coarse parts */
238be712e4SBarry Smith IS coarseparts; /* partitioning on coarse level */
248be712e4SBarry Smith IS fineparts; /* partitioning on fine level */
258be712e4SBarry Smith MatPartitioning coarseMatPart; /* MatPartititioning on coarse level (first level) */
268be712e4SBarry Smith MatPartitioning fineMatPart; /* MatPartitioning on fine level (second level) */
278be712e4SBarry Smith MatPartitioning improver; /* Improve the quality of a partition */
288be712e4SBarry Smith } MatPartitioning_Hierarchical;
298be712e4SBarry Smith
308be712e4SBarry Smith /*
318be712e4SBarry Smith Uses a hierarchical partitioning strategy to partition the matrix in parallel.
328be712e4SBarry Smith Use this interface to make the partitioner consistent with others
338be712e4SBarry Smith */
MatPartitioningApply_Hierarchical(MatPartitioning part,IS * partitioning)348be712e4SBarry Smith static PetscErrorCode MatPartitioningApply_Hierarchical(MatPartitioning part, IS *partitioning)
358be712e4SBarry Smith {
368be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
378be712e4SBarry Smith const PetscInt *fineparts_indices, *coarseparts_indices;
388be712e4SBarry Smith PetscInt *fineparts_indices_tmp;
398be712e4SBarry Smith PetscInt *parts_indices, i, j, mat_localsize, *offsets;
408be712e4SBarry Smith Mat mat = part->adj, adj, sadj;
418be712e4SBarry Smith PetscReal *part_weights;
428be712e4SBarry Smith PetscBool flg;
438be712e4SBarry Smith PetscInt bs = 1;
448be712e4SBarry Smith PetscInt *coarse_vertex_weights = NULL;
458be712e4SBarry Smith PetscMPIInt size, rank;
468be712e4SBarry Smith MPI_Comm comm, scomm;
478be712e4SBarry Smith IS destination, fineparts_temp, vweights, svweights;
488be712e4SBarry Smith PetscInt nsvwegihts, *fp_vweights;
498be712e4SBarry Smith const PetscInt *svweights_indices;
508be712e4SBarry Smith ISLocalToGlobalMapping mapping;
518be712e4SBarry Smith const char *prefix;
528be712e4SBarry Smith PetscBool use_edge_weights;
538be712e4SBarry Smith
548be712e4SBarry Smith PetscFunctionBegin;
558be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
568be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size));
578be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
588be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
598be712e4SBarry Smith if (flg) {
608be712e4SBarry Smith adj = mat;
618be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)adj));
628be712e4SBarry Smith } else {
638be712e4SBarry Smith /* bs indicates if the converted matrix is "reduced" from the original and hence the
648be712e4SBarry Smith resulting partition results need to be stretched to match the original matrix */
658be712e4SBarry Smith PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
668be712e4SBarry Smith if (adj->rmap->n > 0) bs = mat->rmap->n / adj->rmap->n;
678be712e4SBarry Smith }
688be712e4SBarry Smith /* local size of mat */
698be712e4SBarry Smith mat_localsize = adj->rmap->n;
708be712e4SBarry Smith /* check parameters */
718be712e4SBarry Smith /* how many small subdomains we want from a given 'big' suddomain */
728be712e4SBarry Smith PetscCheck(hpart->nfineparts, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, " must set number of small subdomains for each big subdomain ");
738be712e4SBarry 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 ");
748be712e4SBarry Smith
758be712e4SBarry Smith /* Partitioning the domain into one single subdomain is a trivial case, and we should just return */
768be712e4SBarry Smith if (part->n == 1) {
778be712e4SBarry Smith PetscCall(PetscCalloc1(bs * adj->rmap->n, &parts_indices));
788be712e4SBarry Smith PetscCall(ISCreateGeneral(comm, bs * adj->rmap->n, parts_indices, PETSC_OWN_POINTER, partitioning));
798be712e4SBarry Smith hpart->ncoarseparts = 1;
808be712e4SBarry Smith hpart->nfineparts = 1;
818be712e4SBarry Smith PetscCall(PetscStrallocpy("NONE", &hpart->coarseparttype));
828be712e4SBarry Smith PetscCall(PetscStrallocpy("NONE", &hpart->fineparttype));
838be712e4SBarry Smith PetscCall(MatDestroy(&adj));
848be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
858be712e4SBarry Smith }
868be712e4SBarry Smith
878be712e4SBarry Smith if (part->n) {
888be712e4SBarry Smith hpart->ncoarseparts = part->n / hpart->nfineparts;
898be712e4SBarry Smith
908be712e4SBarry Smith if (part->n % hpart->nfineparts != 0) hpart->ncoarseparts++;
918be712e4SBarry Smith } else {
928be712e4SBarry Smith part->n = hpart->ncoarseparts * hpart->nfineparts;
938be712e4SBarry Smith }
948be712e4SBarry Smith
958be712e4SBarry Smith PetscCall(PetscMalloc1(hpart->ncoarseparts + 1, &offsets));
968be712e4SBarry Smith PetscCall(PetscMalloc1(hpart->ncoarseparts, &part_weights));
978be712e4SBarry Smith
988be712e4SBarry Smith offsets[0] = 0;
998be712e4SBarry Smith if (part->n % hpart->nfineparts != 0) offsets[1] = part->n % hpart->nfineparts;
1008be712e4SBarry Smith else offsets[1] = hpart->nfineparts;
1018be712e4SBarry Smith
1028be712e4SBarry Smith part_weights[0] = ((PetscReal)offsets[1]) / part->n;
1038be712e4SBarry Smith
1048be712e4SBarry Smith for (i = 2; i <= hpart->ncoarseparts; i++) {
1058be712e4SBarry Smith offsets[i] = hpart->nfineparts;
1068be712e4SBarry Smith part_weights[i - 1] = ((PetscReal)offsets[i]) / part->n;
1078be712e4SBarry Smith }
1088be712e4SBarry Smith
1098be712e4SBarry Smith offsets[0] = 0;
1108be712e4SBarry Smith for (i = 1; i <= hpart->ncoarseparts; i++) offsets[i] += offsets[i - 1];
1118be712e4SBarry Smith
1128be712e4SBarry Smith /* If these exists a mat partitioner, we should delete it */
1138be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->coarseMatPart));
1148be712e4SBarry Smith PetscCall(MatPartitioningCreate(comm, &hpart->coarseMatPart));
1158be712e4SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)part, &prefix));
1168be712e4SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->coarseMatPart, prefix));
1178be712e4SBarry Smith PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->coarseMatPart, "hierarch_coarse_"));
1188be712e4SBarry Smith /* if did not set partitioning type yet, use parmetis by default */
1198be712e4SBarry Smith if (!hpart->coarseparttype) {
1208be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
1218be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->coarseMatPart, MATPARTITIONINGPARMETIS));
1228be712e4SBarry Smith PetscCall(PetscStrallocpy(MATPARTITIONINGPARMETIS, &hpart->coarseparttype));
1238be712e4SBarry Smith #elif defined(PETSC_HAVE_PTSCOTCH)
1248be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->coarseMatPart, MATPARTITIONINGPTSCOTCH));
1258be712e4SBarry Smith PetscCall(PetscStrallocpy(MATPARTITIONINGPTSCOTCH, &hpart->coarseparttype));
1268be712e4SBarry Smith #else
1278be712e4SBarry Smith SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
1288be712e4SBarry Smith #endif
1298be712e4SBarry Smith } else {
1308be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->coarseMatPart, hpart->coarseparttype));
1318be712e4SBarry Smith }
1328be712e4SBarry Smith PetscCall(MatPartitioningSetAdjacency(hpart->coarseMatPart, adj));
1338be712e4SBarry Smith PetscCall(MatPartitioningSetNParts(hpart->coarseMatPart, hpart->ncoarseparts));
1348be712e4SBarry Smith /* copy over vertex weights */
1358be712e4SBarry Smith if (part->vertex_weights) {
1368be712e4SBarry Smith PetscCall(PetscMalloc1(mat_localsize, &coarse_vertex_weights));
1378be712e4SBarry Smith PetscCall(PetscArraycpy(coarse_vertex_weights, part->vertex_weights, mat_localsize));
1388be712e4SBarry Smith PetscCall(MatPartitioningSetVertexWeights(hpart->coarseMatPart, coarse_vertex_weights));
1398be712e4SBarry Smith }
1408be712e4SBarry Smith /* Copy use_edge_weights flag from part to coarse part */
1418be712e4SBarry Smith PetscCall(MatPartitioningGetUseEdgeWeights(part, &use_edge_weights));
1428be712e4SBarry Smith PetscCall(MatPartitioningSetUseEdgeWeights(hpart->coarseMatPart, use_edge_weights));
1438be712e4SBarry Smith
1448be712e4SBarry Smith PetscCall(MatPartitioningSetPartitionWeights(hpart->coarseMatPart, part_weights));
1458be712e4SBarry Smith PetscCall(MatPartitioningApply(hpart->coarseMatPart, &hpart->coarseparts));
1468be712e4SBarry Smith
1478be712e4SBarry Smith /* Wrap the original vertex weights into an index set so that we can extract the corresponding
1488be712e4SBarry Smith * vertex weights for each big subdomain using ISCreateSubIS().
1498be712e4SBarry Smith * */
1508be712e4SBarry Smith if (part->vertex_weights) PetscCall(ISCreateGeneral(comm, mat_localsize, part->vertex_weights, PETSC_COPY_VALUES, &vweights));
1518be712e4SBarry Smith
1528be712e4SBarry Smith PetscCall(PetscCalloc1(mat_localsize, &fineparts_indices_tmp));
1538be712e4SBarry Smith for (i = 0; i < hpart->ncoarseparts; i += size) {
1548be712e4SBarry Smith /* Determine where we want to send big subdomains */
1558be712e4SBarry Smith PetscCall(MatPartitioningHierarchical_DetermineDestination(part, hpart->coarseparts, i, i + size, &destination));
1568be712e4SBarry Smith /* Assemble a submatrix and its vertex weights for partitioning subdomains */
1578be712e4SBarry Smith PetscCall(MatPartitioningHierarchical_AssembleSubdomain(adj, part->vertex_weights ? vweights : NULL, destination, part->vertex_weights ? &svweights : NULL, &sadj, &mapping));
1588be712e4SBarry Smith /* We have to create a new array to hold vertex weights since coarse partitioner needs to own the vertex-weights array */
1598be712e4SBarry Smith if (part->vertex_weights) {
1608be712e4SBarry Smith PetscCall(ISGetLocalSize(svweights, &nsvwegihts));
1618be712e4SBarry Smith PetscCall(PetscMalloc1(nsvwegihts, &fp_vweights));
1628be712e4SBarry Smith PetscCall(ISGetIndices(svweights, &svweights_indices));
1638be712e4SBarry Smith PetscCall(PetscArraycpy(fp_vweights, svweights_indices, nsvwegihts));
1648be712e4SBarry Smith PetscCall(ISRestoreIndices(svweights, &svweights_indices));
1658be712e4SBarry Smith PetscCall(ISDestroy(&svweights));
1668be712e4SBarry Smith }
1678be712e4SBarry Smith
1688be712e4SBarry Smith PetscCall(ISDestroy(&destination));
1698be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)sadj, &scomm));
1708be712e4SBarry Smith
1718be712e4SBarry Smith /*
1728be712e4SBarry Smith * If the number of big subdomains is smaller than the number of processor cores, the higher ranks do not
1738be712e4SBarry Smith * need to do partitioning
1748be712e4SBarry Smith * */
1758be712e4SBarry Smith if ((i + rank) < hpart->ncoarseparts) {
1768be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->fineMatPart));
1778be712e4SBarry Smith /* create a fine partitioner */
1788be712e4SBarry Smith PetscCall(MatPartitioningCreate(scomm, &hpart->fineMatPart));
1798be712e4SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->fineMatPart, prefix));
1808be712e4SBarry Smith PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->fineMatPart, "hierarch_fine_"));
1818be712e4SBarry Smith /* if do not set partitioning type, use parmetis by default */
1828be712e4SBarry Smith if (!hpart->fineparttype) {
1838be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
1848be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPARMETIS));
1858be712e4SBarry Smith PetscCall(PetscStrallocpy(MATPARTITIONINGPARMETIS, &hpart->fineparttype));
1868be712e4SBarry Smith #elif defined(PETSC_HAVE_PTSCOTCH)
1878be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPTSCOTCH));
1888be712e4SBarry Smith PetscCall(PetscStrallocpy(MATPARTITIONINGPTSCOTCH, &hpart->fineparttype));
1898be712e4SBarry Smith #elif defined(PETSC_HAVE_CHACO)
1908be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGCHACO));
1918be712e4SBarry Smith PetscCall(PetscStrallocpy(MATPARTITIONINGCHACO, &hpart->fineparttype));
1928be712e4SBarry Smith #elif defined(PETSC_HAVE_PARTY)
1938be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->fineMatPart, MATPARTITIONINGPARTY));
1948be712e4SBarry Smith PetscCall(PetscStrallocpy(PETSC_HAVE_PARTY, &hpart->fineparttype));
1958be712e4SBarry Smith #else
1968be712e4SBarry Smith SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis or run with -mat_partitioning_hierarchical_coarseparttype partitiontype");
1978be712e4SBarry Smith #endif
1988be712e4SBarry Smith } else {
1998be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->fineMatPart, hpart->fineparttype));
2008be712e4SBarry Smith }
2018be712e4SBarry Smith PetscCall(MatPartitioningSetUseEdgeWeights(hpart->fineMatPart, use_edge_weights));
2028be712e4SBarry Smith PetscCall(MatPartitioningSetAdjacency(hpart->fineMatPart, sadj));
2038be712e4SBarry Smith PetscCall(MatPartitioningSetNParts(hpart->fineMatPart, offsets[rank + 1 + i] - offsets[rank + i]));
2048be712e4SBarry Smith if (part->vertex_weights) PetscCall(MatPartitioningSetVertexWeights(hpart->fineMatPart, fp_vweights));
2058be712e4SBarry Smith PetscCall(MatPartitioningApply(hpart->fineMatPart, &fineparts_temp));
2068be712e4SBarry Smith } else {
2078be712e4SBarry Smith PetscCall(ISCreateGeneral(scomm, 0, NULL, PETSC_OWN_POINTER, &fineparts_temp));
2088be712e4SBarry Smith }
2098be712e4SBarry Smith
2108be712e4SBarry Smith PetscCall(MatDestroy(&sadj));
2118be712e4SBarry Smith
2128be712e4SBarry Smith /* Send partition back to the original owners */
2138be712e4SBarry Smith PetscCall(MatPartitioningHierarchical_ReassembleFineparts(adj, fineparts_temp, mapping, &hpart->fineparts));
2148be712e4SBarry Smith PetscCall(ISGetIndices(hpart->fineparts, &fineparts_indices));
2158be712e4SBarry Smith for (j = 0; j < mat_localsize; j++)
2168be712e4SBarry Smith if (fineparts_indices[j] >= 0) fineparts_indices_tmp[j] = fineparts_indices[j];
2178be712e4SBarry Smith
2188be712e4SBarry Smith PetscCall(ISRestoreIndices(hpart->fineparts, &fineparts_indices));
2198be712e4SBarry Smith PetscCall(ISDestroy(&hpart->fineparts));
2208be712e4SBarry Smith PetscCall(ISDestroy(&fineparts_temp));
2218be712e4SBarry Smith PetscCall(ISLocalToGlobalMappingDestroy(&mapping));
2228be712e4SBarry Smith }
2238be712e4SBarry Smith
2248be712e4SBarry Smith if (part->vertex_weights) PetscCall(ISDestroy(&vweights));
2258be712e4SBarry Smith
2268be712e4SBarry Smith PetscCall(ISCreateGeneral(comm, mat_localsize, fineparts_indices_tmp, PETSC_OWN_POINTER, &hpart->fineparts));
2278be712e4SBarry Smith PetscCall(ISGetIndices(hpart->fineparts, &fineparts_indices));
2288be712e4SBarry Smith PetscCall(ISGetIndices(hpart->coarseparts, &coarseparts_indices));
2298be712e4SBarry Smith PetscCall(PetscMalloc1(bs * adj->rmap->n, &parts_indices));
2308be712e4SBarry Smith /* Modify the local indices to the global indices by combing the coarse partition and the fine partitions */
2318be712e4SBarry Smith for (i = 0; i < adj->rmap->n; i++) {
2328be712e4SBarry Smith for (j = 0; j < bs; j++) parts_indices[bs * i + j] = fineparts_indices[i] + offsets[coarseparts_indices[i]];
2338be712e4SBarry Smith }
2348be712e4SBarry Smith PetscCall(ISRestoreIndices(hpart->fineparts, &fineparts_indices));
2358be712e4SBarry Smith PetscCall(ISRestoreIndices(hpart->coarseparts, &coarseparts_indices));
2368be712e4SBarry Smith PetscCall(PetscFree(offsets));
2378be712e4SBarry Smith PetscCall(ISCreateGeneral(comm, bs * adj->rmap->n, parts_indices, PETSC_OWN_POINTER, partitioning));
2388be712e4SBarry Smith PetscCall(MatDestroy(&adj));
2398be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
2408be712e4SBarry Smith }
2418be712e4SBarry Smith
MatPartitioningHierarchical_ReassembleFineparts(Mat adj,IS fineparts,ISLocalToGlobalMapping mapping,IS * sfineparts)2428be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_ReassembleFineparts(Mat adj, IS fineparts, ISLocalToGlobalMapping mapping, IS *sfineparts)
2438be712e4SBarry Smith {
2448be712e4SBarry Smith PetscInt *local_indices, *global_indices, *sfineparts_indices, localsize, i;
2458be712e4SBarry Smith const PetscInt *ranges, *fineparts_indices;
2468be712e4SBarry Smith PetscMPIInt rank, *owners;
2478be712e4SBarry Smith MPI_Comm comm;
2488be712e4SBarry Smith PetscLayout rmap;
2498be712e4SBarry Smith PetscSFNode *remote;
2508be712e4SBarry Smith PetscSF sf;
2518be712e4SBarry Smith
2528be712e4SBarry Smith PetscFunctionBegin;
2538be712e4SBarry Smith PetscAssertPointer(sfineparts, 4);
2548be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)adj, &comm));
2558be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
2568be712e4SBarry Smith PetscCall(MatGetLayouts(adj, &rmap, NULL));
2578be712e4SBarry Smith PetscCall(ISGetLocalSize(fineparts, &localsize));
2588be712e4SBarry Smith PetscCall(PetscMalloc2(localsize, &global_indices, localsize, &local_indices));
2598be712e4SBarry Smith for (i = 0; i < localsize; i++) local_indices[i] = i;
2608be712e4SBarry Smith /* map local indices back to global so that we can permulate data globally */
2618be712e4SBarry Smith PetscCall(ISLocalToGlobalMappingApply(mapping, localsize, local_indices, global_indices));
2628be712e4SBarry Smith PetscCall(PetscCalloc1(localsize, &owners));
2638be712e4SBarry Smith /* find owners for global indices */
2648be712e4SBarry Smith for (i = 0; i < localsize; i++) PetscCall(PetscLayoutFindOwner(rmap, global_indices[i], &owners[i]));
2658be712e4SBarry Smith PetscCall(PetscLayoutGetRanges(rmap, &ranges));
2668be712e4SBarry Smith PetscCall(PetscMalloc1(ranges[rank + 1] - ranges[rank], &sfineparts_indices));
2678be712e4SBarry Smith
2688be712e4SBarry Smith for (i = 0; i < ranges[rank + 1] - ranges[rank]; i++) sfineparts_indices[i] = -1;
2698be712e4SBarry Smith
2708be712e4SBarry Smith PetscCall(ISGetIndices(fineparts, &fineparts_indices));
2718be712e4SBarry Smith PetscCall(PetscSFCreate(comm, &sf));
2728be712e4SBarry Smith PetscCall(PetscMalloc1(localsize, &remote));
2738be712e4SBarry Smith for (i = 0; i < localsize; i++) {
2748be712e4SBarry Smith remote[i].rank = owners[i];
2758be712e4SBarry Smith remote[i].index = global_indices[i] - ranges[owners[i]];
2768be712e4SBarry Smith }
2778be712e4SBarry Smith PetscCall(PetscSFSetType(sf, PETSCSFBASIC));
2788be712e4SBarry Smith /* not sure how to add prefix to sf */
2798be712e4SBarry Smith PetscCall(PetscSFSetFromOptions(sf));
2808be712e4SBarry Smith PetscCall(PetscSFSetGraph(sf, localsize, localsize, NULL, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
2818be712e4SBarry Smith PetscCall(PetscSFReduceBegin(sf, MPIU_INT, fineparts_indices, sfineparts_indices, MPI_REPLACE));
2828be712e4SBarry Smith PetscCall(PetscSFReduceEnd(sf, MPIU_INT, fineparts_indices, sfineparts_indices, MPI_REPLACE));
2838be712e4SBarry Smith PetscCall(PetscSFDestroy(&sf));
2848be712e4SBarry Smith PetscCall(ISRestoreIndices(fineparts, &fineparts_indices));
2858be712e4SBarry Smith PetscCall(ISCreateGeneral(comm, ranges[rank + 1] - ranges[rank], sfineparts_indices, PETSC_OWN_POINTER, sfineparts));
2868be712e4SBarry Smith PetscCall(PetscFree2(global_indices, local_indices));
2878be712e4SBarry Smith PetscCall(PetscFree(owners));
2888be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
2898be712e4SBarry Smith }
2908be712e4SBarry Smith
MatPartitioningHierarchical_AssembleSubdomain(Mat adj,IS vweights,IS destination,IS * svweights,Mat * sadj,ISLocalToGlobalMapping * mapping)2918be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_AssembleSubdomain(Mat adj, IS vweights, IS destination, IS *svweights, Mat *sadj, ISLocalToGlobalMapping *mapping)
2928be712e4SBarry Smith {
2938be712e4SBarry Smith IS irows, icols;
2948be712e4SBarry Smith PetscInt irows_ln;
2958be712e4SBarry Smith PetscMPIInt rank;
2968be712e4SBarry Smith const PetscInt *irows_indices;
2978be712e4SBarry Smith MPI_Comm comm;
2988be712e4SBarry Smith
2998be712e4SBarry Smith PetscFunctionBegin;
3008be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)adj, &comm));
3018be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
3028be712e4SBarry Smith /* figure out where data comes from */
3038be712e4SBarry Smith PetscCall(ISBuildTwoSided(destination, NULL, &irows));
3048be712e4SBarry Smith PetscCall(ISDuplicate(irows, &icols));
3058be712e4SBarry Smith PetscCall(ISGetLocalSize(irows, &irows_ln));
3068be712e4SBarry Smith PetscCall(ISGetIndices(irows, &irows_indices));
3078be712e4SBarry Smith PetscCall(ISLocalToGlobalMappingCreate(comm, 1, irows_ln, irows_indices, PETSC_COPY_VALUES, mapping));
3088be712e4SBarry Smith PetscCall(ISRestoreIndices(irows, &irows_indices));
3098be712e4SBarry Smith PetscCall(MatCreateSubMatrices(adj, 1, &irows, &icols, MAT_INITIAL_MATRIX, &sadj));
3108be712e4SBarry Smith if (vweights && svweights) PetscCall(ISCreateSubIS(vweights, irows, svweights));
3118be712e4SBarry Smith PetscCall(ISDestroy(&irows));
3128be712e4SBarry Smith PetscCall(ISDestroy(&icols));
3138be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3148be712e4SBarry Smith }
3158be712e4SBarry Smith
MatPartitioningHierarchical_DetermineDestination(MatPartitioning part,IS partitioning,PetscInt pstart,PetscInt pend,IS * destination)3168be712e4SBarry Smith static PetscErrorCode MatPartitioningHierarchical_DetermineDestination(MatPartitioning part, IS partitioning, PetscInt pstart, PetscInt pend, IS *destination)
3178be712e4SBarry Smith {
3188be712e4SBarry Smith MPI_Comm comm;
3196497c311SBarry Smith PetscMPIInt rank, size;
3206497c311SBarry Smith PetscInt plocalsize, *dest_indices, target;
3218be712e4SBarry Smith const PetscInt *part_indices;
3228be712e4SBarry Smith
3238be712e4SBarry Smith PetscFunctionBegin;
3248be712e4SBarry Smith PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
3258be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank));
3268be712e4SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size));
3278be712e4SBarry 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);
3288be712e4SBarry Smith PetscCheck(pstart <= pend, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, " pstart %" PetscInt_FMT " should be smaller than pend %" PetscInt_FMT, pstart, pend);
3298be712e4SBarry Smith PetscCall(ISGetLocalSize(partitioning, &plocalsize));
3308be712e4SBarry Smith PetscCall(PetscMalloc1(plocalsize, &dest_indices));
3318be712e4SBarry Smith PetscCall(ISGetIndices(partitioning, &part_indices));
3326497c311SBarry Smith for (PetscInt i = 0; i < plocalsize; i++) {
3338be712e4SBarry Smith /* compute target */
3348be712e4SBarry Smith target = part_indices[i] - pstart;
3358be712e4SBarry Smith /* mark out of range entity as -1 */
3368be712e4SBarry Smith if (part_indices[i] < pstart || part_indices[i] >= pend) target = -1;
3378be712e4SBarry Smith dest_indices[i] = target;
3388be712e4SBarry Smith }
3398be712e4SBarry Smith PetscCall(ISCreateGeneral(comm, plocalsize, dest_indices, PETSC_OWN_POINTER, destination));
3408be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3418be712e4SBarry Smith }
3428be712e4SBarry Smith
MatPartitioningView_Hierarchical(MatPartitioning part,PetscViewer viewer)3438be712e4SBarry Smith static PetscErrorCode MatPartitioningView_Hierarchical(MatPartitioning part, PetscViewer viewer)
3448be712e4SBarry Smith {
3458be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
3468be712e4SBarry Smith PetscMPIInt rank;
347*9f196a02SMartin Diehl PetscBool isascii;
3488be712e4SBarry Smith PetscViewer sviewer;
3498be712e4SBarry Smith
3508be712e4SBarry Smith PetscFunctionBegin;
3518be712e4SBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank));
352*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
353*9f196a02SMartin Diehl if (isascii) {
3548be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Number of coarse parts: %" PetscInt_FMT "\n", hpart->ncoarseparts));
3558be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse partitioner: %s\n", hpart->coarseparttype));
3568be712e4SBarry Smith if (hpart->coarseMatPart) {
3578be712e4SBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer));
3588be712e4SBarry Smith PetscCall(MatPartitioningView(hpart->coarseMatPart, viewer));
3598be712e4SBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer));
3608be712e4SBarry Smith }
3618be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Number of fine parts: %" PetscInt_FMT "\n", hpart->nfineparts));
3628be712e4SBarry Smith PetscCall(PetscViewerASCIIPrintf(viewer, " Fine partitioner: %s\n", hpart->fineparttype));
3638be712e4SBarry Smith PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
3648be712e4SBarry Smith if (rank == 0 && hpart->fineMatPart) {
3658be712e4SBarry Smith PetscCall(PetscViewerASCIIPushTab(viewer));
3668be712e4SBarry Smith PetscCall(MatPartitioningView(hpart->fineMatPart, sviewer));
3678be712e4SBarry Smith PetscCall(PetscViewerASCIIPopTab(viewer));
3688be712e4SBarry Smith }
3698be712e4SBarry Smith PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
3708be712e4SBarry Smith }
3718be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3728be712e4SBarry Smith }
3738be712e4SBarry Smith
MatPartitioningHierarchicalGetFineparts(MatPartitioning part,IS * fineparts)3748be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalGetFineparts(MatPartitioning part, IS *fineparts)
3758be712e4SBarry Smith {
3768be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
3778be712e4SBarry Smith
3788be712e4SBarry Smith PetscFunctionBegin;
3798be712e4SBarry Smith *fineparts = hpart->fineparts;
3808be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)hpart->fineparts));
3818be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3828be712e4SBarry Smith }
3838be712e4SBarry Smith
MatPartitioningHierarchicalGetCoarseparts(MatPartitioning part,IS * coarseparts)3848be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalGetCoarseparts(MatPartitioning part, IS *coarseparts)
3858be712e4SBarry Smith {
3868be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
3878be712e4SBarry Smith
3888be712e4SBarry Smith PetscFunctionBegin;
3898be712e4SBarry Smith *coarseparts = hpart->coarseparts;
3908be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)hpart->coarseparts));
3918be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
3928be712e4SBarry Smith }
3938be712e4SBarry Smith
MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning part,PetscInt ncoarseparts)3948be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalSetNcoarseparts(MatPartitioning part, PetscInt ncoarseparts)
3958be712e4SBarry Smith {
3968be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
3978be712e4SBarry Smith
3988be712e4SBarry Smith PetscFunctionBegin;
3998be712e4SBarry Smith hpart->ncoarseparts = ncoarseparts;
4008be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4018be712e4SBarry Smith }
4028be712e4SBarry Smith
MatPartitioningHierarchicalSetNfineparts(MatPartitioning part,PetscInt nfineparts)4038be712e4SBarry Smith PetscErrorCode MatPartitioningHierarchicalSetNfineparts(MatPartitioning part, PetscInt nfineparts)
4048be712e4SBarry Smith {
4058be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
4068be712e4SBarry Smith
4078be712e4SBarry Smith PetscFunctionBegin;
4088be712e4SBarry Smith hpart->nfineparts = nfineparts;
4098be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4108be712e4SBarry Smith }
4118be712e4SBarry Smith
MatPartitioningSetFromOptions_Hierarchical(MatPartitioning part,PetscOptionItems PetscOptionsObject)412ce78bad3SBarry Smith static PetscErrorCode MatPartitioningSetFromOptions_Hierarchical(MatPartitioning part, PetscOptionItems PetscOptionsObject)
4138be712e4SBarry Smith {
4148be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
4158be712e4SBarry Smith char value[1024];
4168be712e4SBarry Smith PetscBool flag = PETSC_FALSE;
4178be712e4SBarry Smith
4188be712e4SBarry Smith PetscFunctionBegin;
4198be712e4SBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Set hierarchical partitioning options");
4208be712e4SBarry Smith PetscCall(PetscOptionsString("-mat_partitioning_hierarchical_coarseparttype", "coarse part type", NULL, NULL, value, sizeof(value), &flag));
4218be712e4SBarry Smith if (flag) {
4228be712e4SBarry Smith PetscCall(PetscFree(hpart->coarseparttype));
4238be712e4SBarry Smith PetscCall(PetscStrallocpy(value, &hpart->coarseparttype));
4248be712e4SBarry Smith }
4258be712e4SBarry Smith PetscCall(PetscOptionsString("-mat_partitioning_hierarchical_fineparttype", "fine part type", NULL, NULL, value, sizeof(value), &flag));
4268be712e4SBarry Smith if (flag) {
4278be712e4SBarry Smith PetscCall(PetscFree(hpart->fineparttype));
4288be712e4SBarry Smith PetscCall(PetscStrallocpy(value, &hpart->fineparttype));
4298be712e4SBarry Smith }
4308be712e4SBarry Smith PetscCall(PetscOptionsInt("-mat_partitioning_hierarchical_ncoarseparts", "number of coarse parts", NULL, hpart->ncoarseparts, &hpart->ncoarseparts, &flag));
4318be712e4SBarry Smith PetscCall(PetscOptionsInt("-mat_partitioning_hierarchical_nfineparts", "number of fine parts", NULL, hpart->nfineparts, &hpart->nfineparts, &flag));
4328be712e4SBarry Smith PetscOptionsHeadEnd();
4338be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4348be712e4SBarry Smith }
4358be712e4SBarry Smith
MatPartitioningDestroy_Hierarchical(MatPartitioning part)4368be712e4SBarry Smith static PetscErrorCode MatPartitioningDestroy_Hierarchical(MatPartitioning part)
4378be712e4SBarry Smith {
4388be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
4398be712e4SBarry Smith
4408be712e4SBarry Smith PetscFunctionBegin;
4418be712e4SBarry Smith if (hpart->coarseparttype) PetscCall(PetscFree(hpart->coarseparttype));
4428be712e4SBarry Smith if (hpart->fineparttype) PetscCall(PetscFree(hpart->fineparttype));
4438be712e4SBarry Smith PetscCall(ISDestroy(&hpart->fineparts));
4448be712e4SBarry Smith PetscCall(ISDestroy(&hpart->coarseparts));
4458be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->coarseMatPart));
4468be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->fineMatPart));
4478be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->improver));
4488be712e4SBarry Smith PetscCall(PetscFree(hpart));
4498be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4508be712e4SBarry Smith }
4518be712e4SBarry Smith
4528be712e4SBarry Smith /*
4538be712e4SBarry Smith Improves the quality of a partition
4548be712e4SBarry Smith */
MatPartitioningImprove_Hierarchical(MatPartitioning part,IS * partitioning)4558be712e4SBarry Smith static PetscErrorCode MatPartitioningImprove_Hierarchical(MatPartitioning part, IS *partitioning)
4568be712e4SBarry Smith {
4578be712e4SBarry Smith MatPartitioning_Hierarchical *hpart = (MatPartitioning_Hierarchical *)part->data;
4588be712e4SBarry Smith Mat mat = part->adj, adj;
4598be712e4SBarry Smith PetscBool flg;
4608be712e4SBarry Smith const char *prefix;
4618be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
4628be712e4SBarry Smith PetscInt *vertex_weights;
4638be712e4SBarry Smith #endif
4648be712e4SBarry Smith
4658be712e4SBarry Smith PetscFunctionBegin;
4668be712e4SBarry Smith PetscCall(PetscObjectTypeCompare((PetscObject)mat, MATMPIADJ, &flg));
4678be712e4SBarry Smith if (flg) {
4688be712e4SBarry Smith adj = mat;
4698be712e4SBarry Smith PetscCall(PetscObjectReference((PetscObject)adj));
4708be712e4SBarry Smith } else {
4718be712e4SBarry Smith /* bs indicates if the converted matrix is "reduced" from the original and hence the
4728be712e4SBarry Smith resulting partition results need to be stretched to match the original matrix */
4738be712e4SBarry Smith PetscCall(MatConvert(mat, MATMPIADJ, MAT_INITIAL_MATRIX, &adj));
4748be712e4SBarry Smith }
4758be712e4SBarry Smith
4768be712e4SBarry Smith /* If there exists a mat partitioner, we should delete it */
4778be712e4SBarry Smith PetscCall(MatPartitioningDestroy(&hpart->improver));
4788be712e4SBarry Smith PetscCall(MatPartitioningCreate(PetscObjectComm((PetscObject)part), &hpart->improver));
4798be712e4SBarry Smith PetscCall(PetscObjectGetOptionsPrefix((PetscObject)part, &prefix));
4808be712e4SBarry Smith PetscCall(PetscObjectSetOptionsPrefix((PetscObject)hpart->improver, prefix));
4818be712e4SBarry Smith PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)hpart->improver, "hierarch_improver_"));
4828be712e4SBarry Smith /* Only parmetis supports to refine a partition */
4838be712e4SBarry Smith #if defined(PETSC_HAVE_PARMETIS)
4848be712e4SBarry Smith PetscCall(MatPartitioningSetType(hpart->improver, MATPARTITIONINGPARMETIS));
4858be712e4SBarry Smith PetscCall(MatPartitioningSetAdjacency(hpart->improver, adj));
4868be712e4SBarry Smith PetscCall(MatPartitioningSetNParts(hpart->improver, part->n));
4878be712e4SBarry Smith /* copy over vertex weights */
4888be712e4SBarry Smith if (part->vertex_weights) {
4898be712e4SBarry Smith PetscCall(PetscMalloc1(adj->rmap->n, &vertex_weights));
4908be712e4SBarry Smith PetscCall(PetscArraycpy(vertex_weights, part->vertex_weights, adj->rmap->n));
4918be712e4SBarry Smith PetscCall(MatPartitioningSetVertexWeights(hpart->improver, vertex_weights));
4928be712e4SBarry Smith }
4938be712e4SBarry Smith PetscCall(MatPartitioningImprove(hpart->improver, partitioning));
4948be712e4SBarry Smith PetscCall(MatDestroy(&adj));
4958be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
4968be712e4SBarry Smith #else
4978be712e4SBarry Smith SETERRQ(PetscObjectComm((PetscObject)adj), PETSC_ERR_SUP, "Requires PETSc be installed with ParMetis");
4988be712e4SBarry Smith #endif
4998be712e4SBarry Smith }
5008be712e4SBarry Smith
5018be712e4SBarry Smith /*MC
5028be712e4SBarry Smith MATPARTITIONINGHIERARCH - Creates a partitioning context via hierarchical partitioning strategy.
5038be712e4SBarry Smith The graph is partitioned into a number of subgraphs, and each subgraph is further split into a few smaller
5048be712e4SBarry Smith subgraphs. The idea can be applied in a recursive manner. It is useful when you want to partition the graph
5058be712e4SBarry Smith into a large number of subgraphs (often more than 10K) since partitions obtained with existing partitioners
5068be712e4SBarry Smith such as ParMETIS and PTScotch are far from ideal. The hierarchical partitioning also tries to avoid off-node
5078be712e4SBarry Smith communication as much as possible for multi-core processor. Another user case for the hierarchical partitioning
508613ce9feSSatish Balay is to improve `PCGASM` convergence by generating multi-process connected subdomain.
5098be712e4SBarry Smith
5108be712e4SBarry Smith Collective
5118be712e4SBarry Smith
5128be712e4SBarry Smith Input Parameter:
5138be712e4SBarry Smith . part - the partitioning context
5148be712e4SBarry Smith
5158be712e4SBarry Smith Options Database Keys:
5168be712e4SBarry Smith + -mat_partitioning_hierarchical_coarseparttype - partitioner type at the first level and parmetis is used by default
5178be712e4SBarry Smith . -mat_partitioning_hierarchical_fineparttype - partitioner type at the second level and parmetis is used by default
5188be712e4SBarry Smith . -mat_partitioning_hierarchical_ncoarseparts - number of subgraphs is required at the first level, which is often the number of compute nodes
5198be712e4SBarry Smith - -mat_partitioning_hierarchical_nfineparts - number of smaller subgraphs for each subgraph, which is often the number of cores per compute node
5208be712e4SBarry Smith
5218be712e4SBarry Smith Level: beginner
5228be712e4SBarry Smith
523613ce9feSSatish Balay Note:
524613ce9feSSatish Balay See {cite}`kong2016highly` and {cite}`kongstognergastonpetersonpermannslaughtermartineau2018`.
5258be712e4SBarry Smith
5268be712e4SBarry Smith .seealso: `MatPartitioningSetType()`, `MatPartitioningType`, `MATPARTITIONINGMETIS`, `MATPARTITIONINGPARMETIS`,
5278be712e4SBarry Smith M*/
5288be712e4SBarry Smith
MatPartitioningCreate_Hierarchical(MatPartitioning part)5298be712e4SBarry Smith PETSC_EXTERN PetscErrorCode MatPartitioningCreate_Hierarchical(MatPartitioning part)
5308be712e4SBarry Smith {
5318be712e4SBarry Smith MatPartitioning_Hierarchical *hpart;
5328be712e4SBarry Smith
5338be712e4SBarry Smith PetscFunctionBegin;
5348be712e4SBarry Smith PetscCall(PetscNew(&hpart));
5358be712e4SBarry Smith part->data = (void *)hpart;
5368be712e4SBarry Smith
5378be712e4SBarry Smith hpart->fineparttype = NULL; /* fine level (second) partitioner */
5388be712e4SBarry Smith hpart->coarseparttype = NULL; /* coarse level (first) partitioner */
5398be712e4SBarry Smith hpart->nfineparts = 1; /* we do not further partition coarse partition any more by default */
5408be712e4SBarry Smith hpart->ncoarseparts = 0; /* number of coarse parts (first level) */
5418be712e4SBarry Smith hpart->coarseparts = NULL;
5428be712e4SBarry Smith hpart->fineparts = NULL;
5438be712e4SBarry Smith hpart->coarseMatPart = NULL;
5448be712e4SBarry Smith hpart->fineMatPart = NULL;
5458be712e4SBarry Smith
5468be712e4SBarry Smith part->ops->apply = MatPartitioningApply_Hierarchical;
5478be712e4SBarry Smith part->ops->view = MatPartitioningView_Hierarchical;
5488be712e4SBarry Smith part->ops->destroy = MatPartitioningDestroy_Hierarchical;
5498be712e4SBarry Smith part->ops->setfromoptions = MatPartitioningSetFromOptions_Hierarchical;
5508be712e4SBarry Smith part->ops->improve = MatPartitioningImprove_Hierarchical;
5518be712e4SBarry Smith PetscFunctionReturn(PETSC_SUCCESS);
5528be712e4SBarry Smith }
553