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