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