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