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