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