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