xref: /petsc/src/dm/impls/plex/tests/ex31.c (revision a69119a591a03a9d906b29c0a4e9802e4d7c9795)
1 static char help[]     = "Test that shared points on interface of partitions can be rebalanced.\n\n";
2 static char FILENAME[] = "ex31.c";
3 
4 #include <petscdmplex.h>
5 #include <petscviewerhdf5.h>
6 #include <petscsf.h>
7 
8 typedef struct {
9   PetscBool parallel;        /* Use ParMetis or Metis */
10   PetscBool useInitialGuess; /* Only active when in parallel, uses RefineKway of ParMetis */
11   PetscInt  entityDepth;     /* depth of the entities to rebalance ( 0 => vertices) */
12 } AppCtx;
13 
14 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
15   PetscFunctionBegin;
16   options->parallel        = PETSC_FALSE;
17   options->useInitialGuess = PETSC_FALSE;
18   options->entityDepth     = 0;
19 
20   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
21   PetscCall(PetscOptionsBoundedInt("-entity_depth", "Depth of the entities to rebalance (0 => vertices)", FILENAME, options->entityDepth, &options->entityDepth, NULL, 0));
22   PetscCall(PetscOptionsBool("-parallel", "Use ParMetis instead of Metis", FILENAME, options->parallel, &options->parallel, NULL));
23   PetscCall(PetscOptionsBool("-use_initial_guess", "Use RefineKway function of ParMetis", FILENAME, options->useInitialGuess, &options->useInitialGuess, NULL));
24   PetscOptionsEnd();
25   PetscFunctionReturn(0);
26 }
27 
28 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm) {
29   PetscFunctionBegin;
30   PetscCall(DMCreate(comm, dm));
31   PetscCall(DMSetType(*dm, DMPLEX));
32   PetscCall(DMSetFromOptions(*dm));
33   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
34   PetscFunctionReturn(0);
35 }
36 
37 int main(int argc, char **argv) {
38   MPI_Comm         comm;
39   DM               dm, dmdist;
40   PetscPartitioner part;
41   AppCtx           user;
42   IS               is = NULL;
43   PetscSection     s = NULL, gsection = NULL;
44   PetscMPIInt      size;
45   PetscSF          sf;
46   PetscInt         pStart, pEnd, p, minBefore, maxBefore, minAfter, maxAfter, gSizeBefore, gSizeAfter;
47   PetscBool        success;
48 
49   PetscFunctionBeginUser;
50   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
51   comm = PETSC_COMM_WORLD;
52   PetscCallMPI(MPI_Comm_size(comm, &size));
53   PetscCall(ProcessOptions(comm, &user));
54   PetscCall(CreateMesh(comm, &user, &dm));
55 
56   /* partition dm using PETSCPARTITIONERPARMETIS */
57   PetscCall(DMPlexGetPartitioner(dm, &part));
58   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)part, "p_"));
59   PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERPARMETIS));
60   PetscCall(PetscPartitionerSetFromOptions(part));
61   PetscCall(PetscSectionCreate(comm, &s));
62   PetscCall(PetscPartitionerDMPlexPartition(part, dm, NULL, s, &is));
63 
64   PetscCall(DMPlexDistribute(dm, 0, NULL, &dmdist));
65   if (dmdist) {
66     PetscCall(DMDestroy(&dm));
67     dm = dmdist;
68   }
69 
70   /* cleanup */
71   PetscCall(PetscSectionDestroy(&s));
72   PetscCall(ISDestroy(&is));
73 
74   /* We make a PetscSection with a DOF on every mesh entity of depth
75    * user.entityDepth, then make a global section and look at its storage size.
76    * We do the same thing after the rebalancing and then assert that the size
77    * remains the same. We also make sure that the balance has improved at least
78    * a little bit compared to the initial decomposition. */
79 
80   if (size > 1) {
81     PetscCall(PetscSectionCreate(comm, &s));
82     PetscCall(PetscSectionSetNumFields(s, 1));
83     PetscCall(PetscSectionSetFieldComponents(s, 0, 1));
84     PetscCall(DMPlexGetDepthStratum(dm, user.entityDepth, &pStart, &pEnd));
85     PetscCall(PetscSectionSetChart(s, pStart, pEnd));
86     for (p = pStart; p < pEnd; ++p) {
87       PetscCall(PetscSectionSetDof(s, p, 1));
88       PetscCall(PetscSectionSetFieldDof(s, p, 0, 1));
89     }
90     PetscCall(PetscSectionSetUp(s));
91     PetscCall(DMGetPointSF(dm, &sf));
92     PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection));
93     PetscCall(PetscSectionGetStorageSize(gsection, &gSizeBefore));
94     minBefore = gSizeBefore;
95     maxBefore = gSizeBefore;
96     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeBefore, 1, MPIU_INT, MPI_SUM, comm));
97     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &minBefore, 1, MPIU_INT, MPI_MIN, comm));
98     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &maxBefore, 1, MPIU_INT, MPI_MAX, comm));
99     PetscCall(PetscSectionDestroy(&gsection));
100   }
101 
102   PetscCall(DMPlexRebalanceSharedPoints(dm, user.entityDepth, user.useInitialGuess, user.parallel, &success));
103 
104   if (size > 1) {
105     PetscCall(PetscSectionCreateGlobalSection(s, sf, PETSC_FALSE, PETSC_FALSE, &gsection));
106     PetscCall(PetscSectionGetStorageSize(gsection, &gSizeAfter));
107     minAfter = gSizeAfter;
108     maxAfter = gSizeAfter;
109     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &gSizeAfter, 1, MPIU_INT, MPI_SUM, comm));
110     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &minAfter, 1, MPIU_INT, MPI_MIN, comm));
111     PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &maxAfter, 1, MPIU_INT, MPI_MAX, comm));
112     PetscCheck(gSizeAfter == gSizeBefore, comm, PETSC_ERR_PLIB, "Global section has not the same size before and after.");
113     PetscCheck(minAfter >= minBefore && maxAfter <= maxBefore && (minAfter > minBefore || maxAfter < maxBefore), comm, PETSC_ERR_PLIB, "DMPlexRebalanceSharedPoints did not improve mesh point balance.");
114     PetscCall(PetscSectionDestroy(&gsection));
115     PetscCall(PetscSectionDestroy(&s));
116   }
117 
118   PetscCall(DMDestroy(&dm));
119   PetscCall(PetscFinalize());
120   return 0;
121 }
122 
123 /*TEST
124 
125   testset:
126     args: -dm_plex_dim 3 -dm_plex_simplex 0
127 
128     test:
129       # rebalance a mesh
130       suffix: 0
131       nsize: {{2 3 4}}
132       requires: parmetis
133       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
134 
135     test:
136       # 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).
137       suffix: 1
138       nsize: {{2 3 4}}
139       requires: parmetis
140       args: -dm_plex_box_faces {{2,3,4  5,4,3  7,11,5}} -entity_depth {{0 1}} -parallel TRUE -use_initial_guess TRUE
141 
142     test:
143       # no-op in serial
144       suffix: 2
145       nsize: {{1}}
146       requires: parmetis
147       args: -dm_plex_box_faces 2,3,4 -entity_depth 0 -parallel FALSE -use_initial_guess FALSE
148 
149 TEST*/
150