1c4762a1bSJed Brown static char help[] = "Test that shared points on interface of partitions can be rebalanced.\n\n";
2c4762a1bSJed Brown static char FILENAME[] = "ex31.c";
3c4762a1bSJed Brown
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown #include <petscviewerhdf5.h>
630602db0SMatthew G. Knepley #include <petscsf.h>
7c4762a1bSJed Brown
8c4762a1bSJed Brown typedef struct {
9c4762a1bSJed Brown PetscBool parallel; /* Use ParMetis or Metis */
10c4762a1bSJed Brown PetscBool useInitialGuess; /* Only active when in parallel, uses RefineKway of ParMetis */
11c4762a1bSJed Brown PetscInt entityDepth; /* depth of the entities to rebalance ( 0 => vertices) */
12c4762a1bSJed Brown } AppCtx;
13c4762a1bSJed Brown
ProcessOptions(MPI_Comm comm,AppCtx * options)14d71ae5a4SJacob Faibussowitsch static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
15d71ae5a4SJacob Faibussowitsch {
16c4762a1bSJed Brown PetscFunctionBegin;
17c4762a1bSJed Brown options->parallel = PETSC_FALSE;
18c4762a1bSJed Brown options->useInitialGuess = PETSC_FALSE;
19c4762a1bSJed Brown options->entityDepth = 0;
2030602db0SMatthew G. Knepley
21d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-entity_depth", "Depth of the entities to rebalance (0 => vertices)", FILENAME, options->entityDepth, &options->entityDepth, NULL, 0));
239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-parallel", "Use ParMetis instead of Metis", FILENAME, options->parallel, &options->parallel, NULL));
249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_initial_guess", "Use RefineKway function of ParMetis", FILENAME, options->useInitialGuess, &options->useInitialGuess, NULL));
25d0609cedSBarry Smith PetscOptionsEnd();
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)29d71ae5a4SJacob Faibussowitsch static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
30d71ae5a4SJacob Faibussowitsch {
31c4762a1bSJed Brown PetscFunctionBegin;
329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm));
339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX));
349566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(*dm));
359566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
37c4762a1bSJed Brown }
38c4762a1bSJed Brown
main(int argc,char ** argv)39d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
40d71ae5a4SJacob Faibussowitsch {
41c4762a1bSJed Brown MPI_Comm comm;
42c4762a1bSJed Brown DM dm, dmdist;
43c4762a1bSJed Brown PetscPartitioner part;
44c4762a1bSJed Brown AppCtx user;
45c4762a1bSJed Brown IS is = NULL;
46c4762a1bSJed Brown PetscSection s = NULL, gsection = NULL;
47c4762a1bSJed Brown PetscMPIInt size;
48c4762a1bSJed Brown PetscSF sf;
49c4762a1bSJed Brown PetscInt pStart, pEnd, p, minBefore, maxBefore, minAfter, maxAfter, gSizeBefore, gSizeAfter;
50c4762a1bSJed Brown PetscBool success;
51c4762a1bSJed Brown
52327415f7SBarry Smith PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
54c4762a1bSJed Brown comm = PETSC_COMM_WORLD;
559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
569566063dSJacob Faibussowitsch PetscCall(ProcessOptions(comm, &user));
579566063dSJacob Faibussowitsch PetscCall(CreateMesh(comm, &user, &dm));
58c4762a1bSJed Brown
59c4762a1bSJed Brown /* partition dm using PETSCPARTITIONERPARMETIS */
609566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part));
619566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "p_"));
629566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS));
639566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part));
649566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &s));
659566063dSJacob Faibussowitsch PetscCall(PetscPartitionerDMPlexPartition(part, dm, NULL, s, &is));
66c4762a1bSJed Brown
679566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist));
68c4762a1bSJed Brown if (dmdist) {
699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
70c4762a1bSJed Brown dm = dmdist;
71c4762a1bSJed Brown }
72c4762a1bSJed Brown
73c4762a1bSJed Brown /* cleanup */
749566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s));
759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* We make a PetscSection with a DOF on every mesh entity of depth
78c4762a1bSJed Brown * user.entityDepth, then make a global section and look at its storage size.
79c4762a1bSJed Brown * We do the same thing after the rebalancing and then assert that the size
80c4762a1bSJed Brown * remains the same. We also make sure that the balance has improved at least
81c4762a1bSJed Brown * a little bit compared to the initial decomposition. */
82c4762a1bSJed Brown
83c4762a1bSJed Brown if (size > 1) {
849566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &s));
859566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(s, 1));
869566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(s, 0, 1));
879566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, user.entityDepth, &pStart, &pEnd));
889566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(s, pStart, pEnd));
89c4762a1bSJed Brown for (p = pStart; p < pEnd; ++p) {
909566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(s, p, 1));
919566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(s, p, 0, 1));
92c4762a1bSJed Brown }
939566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(s));
949566063dSJacob Faibussowitsch PetscCall(DMGetPointSF(dm, &sf));
95eb9d3e4dSMatthew G. Knepley PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &gsection));
969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(gsection, &gSizeBefore));
97c4762a1bSJed Brown minBefore = gSizeBefore;
98c4762a1bSJed Brown maxBefore = gSizeBefore;
99462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gSizeBefore, 1, MPIU_INT, MPI_SUM, comm));
100462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &minBefore, 1, MPIU_INT, MPI_MIN, comm));
101462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &maxBefore, 1, MPIU_INT, MPI_MAX, comm));
1029566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gsection));
103c4762a1bSJed Brown }
104c4762a1bSJed Brown
1059566063dSJacob Faibussowitsch PetscCall(DMPlexRebalanceSharedPoints(dm, user.entityDepth, user.useInitialGuess, user.parallel, &success));
106c4762a1bSJed Brown
107c4762a1bSJed Brown if (size > 1) {
108eb9d3e4dSMatthew G. Knepley PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &gsection));
1099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(gsection, &gSizeAfter));
110c4762a1bSJed Brown minAfter = gSizeAfter;
111c4762a1bSJed Brown maxAfter = gSizeAfter;
112462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &gSizeAfter, 1, MPIU_INT, MPI_SUM, comm));
113462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &minAfter, 1, MPIU_INT, MPI_MIN, comm));
114462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &maxAfter, 1, MPIU_INT, MPI_MAX, comm));
11508401ef6SPierre Jolivet PetscCheck(gSizeAfter == gSizeBefore, comm, PETSC_ERR_PLIB, "Global section has not the same size before and after.");
1161dca8a05SBarry Smith PetscCheck(minAfter >= minBefore && maxAfter <= maxBefore && (minAfter > minBefore || maxAfter < maxBefore), comm, PETSC_ERR_PLIB, "DMPlexRebalanceSharedPoints did not improve mesh point balance.");
1179566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&gsection));
1189566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&s));
119c4762a1bSJed Brown }
120c4762a1bSJed Brown
1219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm));
1229566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
123b122ec5aSJacob Faibussowitsch return 0;
124c4762a1bSJed Brown }
125c4762a1bSJed Brown
126c4762a1bSJed Brown /*TEST
127c4762a1bSJed Brown
12830602db0SMatthew G. Knepley testset:
12930602db0SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0
13030602db0SMatthew G. Knepley
131c4762a1bSJed Brown test:
132c4762a1bSJed Brown # rebalance a mesh
133c4762a1bSJed Brown suffix: 0
134c4762a1bSJed Brown nsize: {{2 3 4}}
135c4762a1bSJed Brown requires: parmetis
13630602db0SMatthew G. Knepley args: -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -entity_depth {{0 1}} -parallel {{FALSE TRUE}} -use_initial_guess FALSE
137*e0008caeSPierre Jolivet output_file: output/empty.out
138c4762a1bSJed Brown
139c4762a1bSJed Brown test:
140c4762a1bSJed Brown # rebalance a mesh but use the initial guess (uses a random algorithm and gives different results on different machines, so just check that it runs).
141c4762a1bSJed Brown suffix: 1
142c4762a1bSJed Brown nsize: {{2 3 4}}
143c4762a1bSJed Brown requires: parmetis
14430602db0SMatthew G. Knepley args: -dm_plex_box_faces {{2,3,4 5,4,3 7,11,5}} -entity_depth {{0 1}} -parallel TRUE -use_initial_guess TRUE
1453886731fSPierre Jolivet output_file: output/empty.out
146c4762a1bSJed Brown
147c4762a1bSJed Brown test:
148c4762a1bSJed Brown # no-op in serial
149c4762a1bSJed Brown suffix: 2
150c4762a1bSJed Brown nsize: {{1}}
151c4762a1bSJed Brown requires: parmetis
15230602db0SMatthew G. Knepley args: -dm_plex_box_faces 2,3,4 -entity_depth 0 -parallel FALSE -use_initial_guess FALSE
1533886731fSPierre Jolivet output_file: output/empty.out
154c4762a1bSJed Brown
155c4762a1bSJed Brown TEST*/
156