xref: /petsc/src/dm/impls/plex/tests/ex24.c (revision bef158480efac06de457f7a665168877ab3c2fd7) !
1 static char help[] = "Test that MatPartitioning and PetscPartitioner interfaces are equivalent when using PETSCPARTITIONERMATPARTITIONING\n\n";
2 static char FILENAME[] = "ex24.c";
3 
4 #include <petscdmplex.h>
5 #include <petscviewerhdf5.h>
6 
7 #if defined(PETSC_HAVE_PTSCOTCH)
8 EXTERN_C_BEGIN
9 #include <ptscotch.h>
10 EXTERN_C_END
11 #endif
12 
13 typedef struct {
14   PetscInt  dim;                          /* The topological mesh dimension */
15   PetscInt  faces[3];                     /* Number of faces per dimension */
16   PetscBool simplex;                      /* Use simplices or hexes */
17   PetscBool interpolate;                  /* Interpolate mesh */
18   PetscBool compare_is;                   /* Compare ISs and PetscSections */
19   PetscBool compare_dm;                   /* Compare DM */
20   PetscBool tpw;                          /* Use target partition weights */
21   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
22   char      partitioning[64];
23   char      repartitioning[64];
24 } AppCtx;
25 
26 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
27 {
28   PetscInt dim;
29   PetscBool repartition = PETSC_TRUE;
30   PetscErrorCode ierr;
31 
32   PetscFunctionBegin;
33   options->compare_is   = PETSC_FALSE;
34   options->compare_dm   = PETSC_FALSE;
35   options->dim          = 3;
36   options->simplex      = PETSC_TRUE;
37   options->interpolate  = PETSC_FALSE;
38   options->filename[0]  = '\0';
39   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
40   ierr = PetscOptionsBool("-compare_is", "Compare ISs and PetscSections?", FILENAME, options->compare_is, &options->compare_is, NULL);CHKERRQ(ierr);
41   ierr = PetscOptionsBool("-compare_dm", "Compare DMs?", FILENAME, options->compare_dm, &options->compare_dm, NULL);CHKERRQ(ierr);
42   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", FILENAME, options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
43   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", FILENAME, options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
44   ierr = PetscOptionsBool("-interpolate", "Interpolate the mesh", FILENAME, options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
45   ierr = PetscOptionsString("-filename", "The mesh file", FILENAME, options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
46   options->faces[0] = 1; options->faces[1] = 1; options->faces[2] = 1;
47   dim = options->dim;
48   ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", FILENAME, options->faces, &dim, NULL);CHKERRQ(ierr);
49   if (dim) options->dim = dim;
50   ierr = PetscStrncpy(options->partitioning,MATPARTITIONINGPARMETIS,sizeof(options->partitioning));CHKERRQ(ierr);
51   ierr = PetscOptionsString("-partitioning","The mat partitioning type to test","None",options->partitioning, options->partitioning,sizeof(options->partitioning),NULL);CHKERRQ(ierr);
52   ierr = PetscOptionsBool("-repartition", "Partition again after the first partition?", FILENAME, repartition, &repartition, NULL);CHKERRQ(ierr);
53   if (repartition) {
54     ierr = PetscStrncpy(options->repartitioning,MATPARTITIONINGPARMETIS,64);CHKERRQ(ierr);
55     ierr = PetscOptionsString("-repartitioning","The mat partitioning type to test (second partitioning)","None", options->repartitioning, options->repartitioning,sizeof(options->repartitioning),NULL);CHKERRQ(ierr);
56   } else {
57     options->repartitioning[0] = '\0';
58   }
59   ierr = PetscOptionsBool("-tpweight", "Use target partition weights", FILENAME, options->tpw, &options->tpw, NULL);CHKERRQ(ierr);
60   ierr = PetscOptionsEnd();
61   PetscFunctionReturn(0);
62 }
63 
64 static PetscErrorCode ScotchResetRandomSeed()
65 {
66 #if defined(PETSC_HAVE_PTSCOTCH)
67   SCOTCH_randomReset();
68 #endif
69   PetscFunctionReturn(0);
70 }
71 
72 
73 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
74 {
75   PetscInt       dim          = user->dim;
76   PetscInt      *faces        = user->faces;
77   PetscBool      simplex      = user->simplex;
78   PetscBool      interpolate  = user->interpolate;
79   const char    *filename     = user->filename;
80   size_t         len;
81   PetscMPIInt    rank;
82   PetscErrorCode ierr;
83 
84   PetscFunctionBegin;
85   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
86   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
87   if (len) {
88     ierr = DMPlexCreateFromFile(comm, filename, interpolate, dm);CHKERRQ(ierr);
89     ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
90   } else {
91     ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, interpolate, dm);CHKERRQ(ierr);
92   }
93   PetscFunctionReturn(0);
94 }
95 
96 int main(int argc, char **argv)
97 {
98   MPI_Comm       comm;
99   DM             dm1, dm2, dmdist1, dmdist2;
100   MatPartitioning mp;
101   PetscPartitioner part1, part2;
102   AppCtx         user;
103   IS             is1=NULL, is2=NULL;
104   IS             is1g, is2g;
105   PetscSection   s1=NULL, s2=NULL, tpws = NULL;
106   PetscInt       i;
107   PetscBool      flg;
108   PetscErrorCode ierr;
109   PetscMPIInt    size;
110 
111   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
112   comm = PETSC_COMM_WORLD;
113   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
114   ierr = ProcessOptions(comm, &user);CHKERRQ(ierr);
115   ierr = CreateMesh(comm, &user, &dm1);CHKERRQ(ierr);
116   ierr = CreateMesh(comm, &user, &dm2);CHKERRQ(ierr);
117 
118   if (user.tpw) {
119     ierr = PetscSectionCreate(comm, &tpws);CHKERRQ(ierr);
120     ierr = PetscSectionSetChart(tpws, 0, size);CHKERRQ(ierr);
121     for (i=0;i<size;i++) {
122       PetscInt tdof = i%2 ? 2*i -1 : i+2;
123       ierr = PetscSectionSetDof(tpws, i, tdof);CHKERRQ(ierr);
124     }
125     if (size > 1) { /* test zero tpw entry */
126       ierr = PetscSectionSetDof(tpws, 0, 0);CHKERRQ(ierr);
127     }
128     ierr = PetscSectionSetUp(tpws);CHKERRQ(ierr);
129   }
130 
131   /* partition dm1 using PETSCPARTITIONERPARMETIS */
132   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
133   ierr = DMPlexGetPartitioner(dm1, &part1);CHKERRQ(ierr);
134   ierr = PetscObjectSetOptionsPrefix((PetscObject)part1,"p1_");CHKERRQ(ierr);
135   ierr = PetscPartitionerSetType(part1, user.partitioning);CHKERRQ(ierr);
136   ierr = PetscPartitionerSetFromOptions(part1);CHKERRQ(ierr);
137   ierr = PetscSectionCreate(comm, &s1);CHKERRQ(ierr);
138   ierr = PetscPartitionerDMPlexPartition(part1, dm1, tpws, s1, &is1);CHKERRQ(ierr);
139 
140   /* partition dm2 using PETSCPARTITIONERMATPARTITIONING with MATPARTITIONINGPARMETIS */
141   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
142   ierr = DMPlexGetPartitioner(dm2, &part2);CHKERRQ(ierr);
143   ierr = PetscObjectSetOptionsPrefix((PetscObject)part2,"p2_");CHKERRQ(ierr);
144   ierr = PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING);CHKERRQ(ierr);
145   ierr = PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp);CHKERRQ(ierr);
146   ierr = MatPartitioningSetType(mp, user.partitioning);CHKERRQ(ierr);
147   ierr = PetscPartitionerSetFromOptions(part2);CHKERRQ(ierr);
148   ierr = PetscSectionCreate(comm, &s2);CHKERRQ(ierr);
149   ierr = PetscPartitionerDMPlexPartition(part2, dm2, tpws, s2, &is2);CHKERRQ(ierr);
150 
151   ierr = ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g);CHKERRQ(ierr);
152   ierr = ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g);CHKERRQ(ierr);
153   ierr = ISViewFromOptions(is1g, NULL, "-seq_is1_view");CHKERRQ(ierr);
154   ierr = ISViewFromOptions(is2g, NULL, "-seq_is2_view");CHKERRQ(ierr);
155   /* compare the two ISs */
156   if (user.compare_is) {
157     ierr = ISEqualUnsorted(is1g, is2g, &flg);CHKERRQ(ierr);
158     if (!flg) {ierr = PetscPrintf(comm, "ISs are not equal with type %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);}
159   }
160   ierr = ISDestroy(&is1g);CHKERRQ(ierr);
161   ierr = ISDestroy(&is2g);CHKERRQ(ierr);
162 
163   /* compare the two PetscSections */
164   ierr = PetscSectionViewFromOptions(s1, NULL, "-seq_s1_view");CHKERRQ(ierr);
165   ierr = PetscSectionViewFromOptions(s2, NULL, "-seq_s2_view");CHKERRQ(ierr);
166   if (user.compare_is) {
167     ierr = PetscSectionCompare(s1, s2, &flg);CHKERRQ(ierr);
168     if (!flg) {ierr = PetscPrintf(comm, "PetscSections are not equal with %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);}
169   }
170 
171   /* distribute both DMs */
172   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
173   ierr = DMPlexDistribute(dm1, 0, NULL, &dmdist1);CHKERRQ(ierr);
174   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
175   ierr = DMPlexDistribute(dm2, 0, NULL, &dmdist2);CHKERRQ(ierr);
176 
177   /* cleanup */
178   ierr = PetscSectionDestroy(&tpws);CHKERRQ(ierr);
179   ierr = PetscSectionDestroy(&s1);CHKERRQ(ierr);
180   ierr = PetscSectionDestroy(&s2);CHKERRQ(ierr);
181   ierr = ISDestroy(&is1);CHKERRQ(ierr);
182   ierr = ISDestroy(&is2);CHKERRQ(ierr);
183   ierr = DMDestroy(&dm1);CHKERRQ(ierr);
184   ierr = DMDestroy(&dm2);CHKERRQ(ierr);
185 
186   /* if distributed DMs are NULL (sequential case), then quit */
187   if (!dmdist1 && !dmdist2) return ierr;
188 
189   ierr = DMViewFromOptions(dmdist1, NULL, "-dm_dist1_view");CHKERRQ(ierr);
190   ierr = DMViewFromOptions(dmdist2, NULL, "-dm_dist2_view");CHKERRQ(ierr);
191 
192   /* compare the two distributed DMs */
193   if (user.compare_dm) {
194     ierr = DMPlexEqual(dmdist1, dmdist2, &flg);CHKERRQ(ierr);
195     if (!flg) {ierr = PetscPrintf(comm, "Distributed DMs are not equal %s with size %d.\n",user.partitioning,size);CHKERRQ(ierr);}
196   }
197 
198   /* if repartitioning is disabled, then quit */
199   if (user.repartitioning[0] == '\0') return ierr;
200 
201   if (user.tpw) {
202     ierr = PetscSectionCreate(comm, &tpws);CHKERRQ(ierr);
203     ierr = PetscSectionSetChart(tpws, 0, size);CHKERRQ(ierr);
204     for (i=0;i<size;i++) {
205       PetscInt tdof = i%2 ? i+1 : size - i;
206       ierr = PetscSectionSetDof(tpws, i, tdof);CHKERRQ(ierr);
207     }
208     ierr = PetscSectionSetUp(tpws);CHKERRQ(ierr);
209   }
210 
211   /* repartition distributed DM dmdist1 */
212   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
213   ierr = DMPlexGetPartitioner(dmdist1, &part1);CHKERRQ(ierr);
214   ierr = PetscObjectSetOptionsPrefix((PetscObject)part1,"dp1_");CHKERRQ(ierr);
215   ierr = PetscPartitionerSetType(part1, user.repartitioning);CHKERRQ(ierr);
216   ierr = PetscPartitionerSetFromOptions(part1);CHKERRQ(ierr);
217   ierr = PetscSectionCreate(comm, &s1);CHKERRQ(ierr);
218   ierr = PetscPartitionerDMPlexPartition(part1, dmdist1, tpws, s1, &is1);CHKERRQ(ierr);
219 
220   /* repartition distributed DM dmdist2 */
221   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
222   ierr = DMPlexGetPartitioner(dmdist2, &part2);CHKERRQ(ierr);
223   ierr = PetscObjectSetOptionsPrefix((PetscObject)part2,"dp2_");CHKERRQ(ierr);
224   ierr = PetscPartitionerSetType(part2, PETSCPARTITIONERMATPARTITIONING);CHKERRQ(ierr);
225   ierr = PetscPartitionerMatPartitioningGetMatPartitioning(part2, &mp);CHKERRQ(ierr);
226   ierr = MatPartitioningSetType(mp, user.repartitioning);CHKERRQ(ierr);
227   ierr = PetscPartitionerSetFromOptions(part2);CHKERRQ(ierr);
228   ierr = PetscSectionCreate(comm, &s2);CHKERRQ(ierr);
229   ierr = PetscPartitionerDMPlexPartition(part2, dmdist2, tpws, s2, &is2);CHKERRQ(ierr);
230 
231   /* compare the two ISs */
232   ierr = ISOnComm(is1, comm, PETSC_USE_POINTER, &is1g);CHKERRQ(ierr);
233   ierr = ISOnComm(is2, comm, PETSC_USE_POINTER, &is2g);CHKERRQ(ierr);
234   ierr = ISViewFromOptions(is1g, NULL, "-dist_is1_view");CHKERRQ(ierr);
235   ierr = ISViewFromOptions(is2g, NULL, "-dist_is2_view");CHKERRQ(ierr);
236   if (user.compare_is) {
237     ierr = ISEqualUnsorted(is1g, is2g, &flg);CHKERRQ(ierr);
238     if (!flg) {ierr = PetscPrintf(comm, "Distributed ISs are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);}
239   }
240   ierr = ISDestroy(&is1g);CHKERRQ(ierr);
241   ierr = ISDestroy(&is2g);CHKERRQ(ierr);
242 
243   /* compare the two PetscSections */
244   ierr = PetscSectionViewFromOptions(s1, NULL, "-dist_s1_view");CHKERRQ(ierr);
245   ierr = PetscSectionViewFromOptions(s2, NULL, "-dist_s2_view");CHKERRQ(ierr);
246   if (user.compare_is) {
247     ierr = PetscSectionCompare(s1, s2, &flg);CHKERRQ(ierr);
248     if (!flg) {ierr = PetscPrintf(comm, "Distributed PetscSections are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);}
249   }
250 
251   /* redistribute both distributed DMs */
252   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
253   ierr = DMPlexDistribute(dmdist1, 0, NULL, &dm1);CHKERRQ(ierr);
254   ierr = ScotchResetRandomSeed();CHKERRQ(ierr);
255   ierr = DMPlexDistribute(dmdist2, 0, NULL, &dm2);CHKERRQ(ierr);
256 
257   /* compare the two distributed DMs */
258   if (!user.interpolate) {
259     ierr = DMPlexEqual(dm1, dm2, &flg);CHKERRQ(ierr);
260     if (!flg) {ierr = PetscPrintf(comm, "Redistributed DMs are not equal, with %s with size %d.\n",user.repartitioning,size);CHKERRQ(ierr);}
261   }
262 
263   /* cleanup */
264   ierr = PetscSectionDestroy(&tpws);CHKERRQ(ierr);
265   ierr = PetscSectionDestroy(&s1);CHKERRQ(ierr);
266   ierr = PetscSectionDestroy(&s2);CHKERRQ(ierr);
267   ierr = ISDestroy(&is1);CHKERRQ(ierr);
268   ierr = ISDestroy(&is2);CHKERRQ(ierr);
269   ierr = DMDestroy(&dm1);CHKERRQ(ierr);
270   ierr = DMDestroy(&dm2);CHKERRQ(ierr);
271   ierr = DMDestroy(&dmdist1);CHKERRQ(ierr);
272   ierr = DMDestroy(&dmdist2);CHKERRQ(ierr);
273   ierr = PetscFinalize();
274   return ierr;
275 }
276 
277 /*TEST
278 
279   test:
280     # partition sequential mesh loaded from Exodus file
281     suffix: 0
282     nsize: {{1 2 3 4 8}}
283     requires: chaco parmetis ptscotch exodusii
284     args: -filename ${PETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -interpolate
285     args: -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -tpweight {{0 1}}
286   test:
287     # repartition mesh already partitioned naively by MED loader
288     suffix: 1
289     nsize: {{1 2 3 4 8}}
290     TODO: MED
291     requires: parmetis ptscotch med
292     args: -filename ${PETSC_DIR}/share/petsc/datafiles/meshes/cylinder.med -interpolate
293     args: -repartition 0 -partitioning {{parmetis ptscotch}}
294   test:
295     # partition mesh generated by ctetgen using scotch, then repartition with scotch, diff view
296     suffix: 3
297     nsize: 4
298     requires: ptscotch ctetgen
299     args: -faces 2,3,2 -partitioning ptscotch -repartitioning ptscotch -interpolate
300     args: -p1_petscpartitioner_view -p2_petscpartitioner_view -dp1_petscpartitioner_view -dp2_petscpartitioner_view -tpweight {{0 1}}
301   test:
302     # partition mesh generated by ctetgen using partitioners supported both by MatPartitioning and PetscPartitioner
303     suffix: 4
304     nsize: {{1 2 3 4 8}}
305     requires: chaco parmetis ptscotch ctetgen
306     args: -faces {{2,3,4  5,4,3  7,11,5}} -partitioning {{chaco parmetis ptscotch}} -repartitioning {{parmetis ptscotch}} -interpolate -tpweight {{0 1}}
307 
308 TEST*/
309 
310