xref: /petsc/src/dm/impls/plex/tests/ex18.c (revision 1dca8a0504492127e77eac64bc165d7372dd6d63)
1c4762a1bSJed Brown static char help[] = "Tests for parallel mesh loading and parallel topological interpolation\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>
4c4762a1bSJed Brown /* List of test meshes
5c4762a1bSJed Brown 
6c4762a1bSJed Brown Network
7c4762a1bSJed Brown -------
8c4762a1bSJed Brown Test 0 (2 ranks):
9c4762a1bSJed Brown 
10c4762a1bSJed Brown network=0:
11c4762a1bSJed Brown ---------
12c4762a1bSJed Brown   cell 0   cell 1   cell 2          nCells-1       (edge)
13c4762a1bSJed Brown 0 ------ 1 ------ 2 ------ 3 -- -- v --  -- nCells (vertex)
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   vertex distribution:
16c4762a1bSJed Brown     rank 0: 0 1
17c4762a1bSJed Brown     rank 1: 2 3 ... nCells
18c4762a1bSJed Brown   cell(edge) distribution:
19c4762a1bSJed Brown     rank 0: 0 1
20c4762a1bSJed Brown     rank 1: 2 ... nCells-1
21c4762a1bSJed Brown 
22c4762a1bSJed Brown network=1:
23c4762a1bSJed Brown ---------
24c4762a1bSJed Brown                v2
25c4762a1bSJed Brown                 ^
26c4762a1bSJed Brown                 |
27c4762a1bSJed Brown                cell 2
28c4762a1bSJed Brown                 |
29c4762a1bSJed Brown  v0 --cell 0--> v3--cell 1--> v1
30c4762a1bSJed Brown 
31c4762a1bSJed Brown   vertex distribution:
32c4762a1bSJed Brown     rank 0: 0 1 3
33c4762a1bSJed Brown     rank 1: 2
34c4762a1bSJed Brown   cell(edge) distribution:
35c4762a1bSJed Brown     rank 0: 0 1
36c4762a1bSJed Brown     rank 1: 2
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   example:
39c4762a1bSJed Brown     mpiexec -n 2 ./ex18 -distribute 1 -dim 1 -orig_dm_view -dist_dm_view -dist_dm_view -petscpartitioner_type parmetis -ncells 50
40c4762a1bSJed Brown 
41c4762a1bSJed Brown Triangle
42c4762a1bSJed Brown --------
43c4762a1bSJed Brown Test 0 (2 ranks):
44c4762a1bSJed Brown Two triangles sharing a face
45c4762a1bSJed Brown 
46c4762a1bSJed Brown         2
47c4762a1bSJed Brown       / | \
48c4762a1bSJed Brown      /  |  \
49c4762a1bSJed Brown     /   |   \
50c4762a1bSJed Brown    0  0 | 1  3
51c4762a1bSJed Brown     \   |   /
52c4762a1bSJed Brown      \  |  /
53c4762a1bSJed Brown       \ | /
54c4762a1bSJed Brown         1
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   vertex distribution:
57c4762a1bSJed Brown     rank 0: 0 1
58c4762a1bSJed Brown     rank 1: 2 3
59c4762a1bSJed Brown   cell distribution:
60c4762a1bSJed Brown     rank 0: 0
61c4762a1bSJed Brown     rank 1: 1
62c4762a1bSJed Brown 
63c4762a1bSJed Brown Test 1 (3 ranks):
64c4762a1bSJed Brown Four triangles partitioned across 3 ranks
65c4762a1bSJed Brown 
66c4762a1bSJed Brown    0 _______ 3
67c4762a1bSJed Brown    | \     / |
68c4762a1bSJed Brown    |  \ 1 /  |
69c4762a1bSJed Brown    |   \ /   |
70c4762a1bSJed Brown    | 0  2  2 |
71c4762a1bSJed Brown    |   / \   |
72c4762a1bSJed Brown    |  / 3 \  |
73c4762a1bSJed Brown    | /     \ |
74c4762a1bSJed Brown    1 ------- 4
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   vertex distribution:
77c4762a1bSJed Brown     rank 0: 0 1
78c4762a1bSJed Brown     rank 1: 2 3
79c4762a1bSJed Brown     rank 2: 4
80c4762a1bSJed Brown   cell distribution:
81c4762a1bSJed Brown     rank 0: 0
82c4762a1bSJed Brown     rank 1: 1
83c4762a1bSJed Brown     rank 2: 2 3
84c4762a1bSJed Brown 
85c4762a1bSJed Brown Test 2 (3 ranks):
86c4762a1bSJed Brown Four triangles partitioned across 3 ranks
87c4762a1bSJed Brown 
88c4762a1bSJed Brown    1 _______ 3
89c4762a1bSJed Brown    | \     / |
90c4762a1bSJed Brown    |  \ 1 /  |
91c4762a1bSJed Brown    |   \ /   |
92c4762a1bSJed Brown    | 0  0  2 |
93c4762a1bSJed Brown    |   / \   |
94c4762a1bSJed Brown    |  / 3 \  |
95c4762a1bSJed Brown    | /     \ |
96c4762a1bSJed Brown    2 ------- 4
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   vertex distribution:
99c4762a1bSJed Brown     rank 0: 0 1
100c4762a1bSJed Brown     rank 1: 2 3
101c4762a1bSJed Brown     rank 2: 4
102c4762a1bSJed Brown   cell distribution:
103c4762a1bSJed Brown     rank 0: 0
104c4762a1bSJed Brown     rank 1: 1
105c4762a1bSJed Brown     rank 2: 2 3
106c4762a1bSJed Brown 
107c4762a1bSJed Brown Tetrahedron
108c4762a1bSJed Brown -----------
109c4762a1bSJed Brown Test 0:
110c4762a1bSJed Brown Two tets sharing a face
111c4762a1bSJed Brown 
112c4762a1bSJed Brown  cell   3 _______    cell
113c4762a1bSJed Brown  0    / | \      \   1
114c4762a1bSJed Brown      /  |  \      \
115c4762a1bSJed Brown     /   |   \      \
116c4762a1bSJed Brown    0----|----4-----2
117c4762a1bSJed Brown     \   |   /      /
118c4762a1bSJed Brown      \  |  /      /
119c4762a1bSJed Brown       \ | /      /
120c4762a1bSJed Brown         1-------
121c4762a1bSJed Brown    y
122c4762a1bSJed Brown    | x
123c4762a1bSJed Brown    |/
124c4762a1bSJed Brown    *----z
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   vertex distribution:
127c4762a1bSJed Brown     rank 0: 0 1
128c4762a1bSJed Brown     rank 1: 2 3 4
129c4762a1bSJed Brown   cell distribution:
130c4762a1bSJed Brown     rank 0: 0
131c4762a1bSJed Brown     rank 1: 1
132c4762a1bSJed Brown 
133c4762a1bSJed Brown Quadrilateral
134c4762a1bSJed Brown -------------
135c4762a1bSJed Brown Test 0 (2 ranks):
136c4762a1bSJed Brown Two quads sharing a face
137c4762a1bSJed Brown 
138c4762a1bSJed Brown    3-------2-------5
139c4762a1bSJed Brown    |       |       |
140c4762a1bSJed Brown    |   0   |   1   |
141c4762a1bSJed Brown    |       |       |
142c4762a1bSJed Brown    0-------1-------4
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   vertex distribution:
145c4762a1bSJed Brown     rank 0: 0 1 2
146c4762a1bSJed Brown     rank 1: 3 4 5
147c4762a1bSJed Brown   cell distribution:
148c4762a1bSJed Brown     rank 0: 0
149c4762a1bSJed Brown     rank 1: 1
150c4762a1bSJed Brown 
151c4762a1bSJed Brown TODO Test 1:
152c4762a1bSJed Brown A quad and a triangle sharing a face
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    5-------4
155c4762a1bSJed Brown    |       | \
156c4762a1bSJed Brown    |   0   |  \
157c4762a1bSJed Brown    |       | 1 \
158c4762a1bSJed Brown    2-------3----6
159c4762a1bSJed Brown 
160c4762a1bSJed Brown Hexahedron
161c4762a1bSJed Brown ----------
162c4762a1bSJed Brown Test 0 (2 ranks):
163c4762a1bSJed Brown Two hexes sharing a face
164c4762a1bSJed Brown 
165c4762a1bSJed Brown cell   7-------------6-------------11 cell
166c4762a1bSJed Brown 0     /|            /|            /|     1
167c4762a1bSJed Brown      / |   F1      / |   F7      / |
168c4762a1bSJed Brown     /  |          /  |          /  |
169c4762a1bSJed Brown    4-------------5-------------10  |
170c4762a1bSJed Brown    |   |     F4  |   |     F10 |   |
171c4762a1bSJed Brown    |   |         |   |         |   |
172c4762a1bSJed Brown    |F5 |         |F3 |         |F9 |
173c4762a1bSJed Brown    |   |  F2     |   |   F8    |   |
174c4762a1bSJed Brown    |   3---------|---2---------|---9
175c4762a1bSJed Brown    |  /          |  /          |  /
176c4762a1bSJed Brown    | /   F0      | /    F6     | /
177c4762a1bSJed Brown    |/            |/            |/
178c4762a1bSJed Brown    0-------------1-------------8
179c4762a1bSJed Brown 
180c4762a1bSJed Brown   vertex distribution:
181c4762a1bSJed Brown     rank 0: 0 1 2 3 4 5
182c4762a1bSJed Brown     rank 1: 6 7 8 9 10 11
183c4762a1bSJed Brown   cell distribution:
184c4762a1bSJed Brown     rank 0: 0
185c4762a1bSJed Brown     rank 1: 1
186c4762a1bSJed Brown 
187c4762a1bSJed Brown */
188c4762a1bSJed Brown 
189c4762a1bSJed Brown typedef enum {NONE, CREATE, AFTER_CREATE, AFTER_DISTRIBUTE} InterpType;
190c4762a1bSJed Brown 
191c4762a1bSJed Brown typedef struct {
192c4762a1bSJed Brown   PetscInt   debug;                        /* The debugging level */
193c4762a1bSJed Brown   PetscInt   testNum;                      /* Indicates the mesh to create */
194c4762a1bSJed Brown   PetscInt   dim;                          /* The topological mesh dimension */
195c4762a1bSJed Brown   PetscBool  cellSimplex;                  /* Use simplices or hexes */
196c4762a1bSJed Brown   PetscBool  distribute;                   /* Distribute the mesh */
197c4762a1bSJed Brown   InterpType interpolate;                  /* Interpolate the mesh before or after DMPlexDistribute() */
198c4762a1bSJed Brown   PetscBool  useGenerator;                 /* Construct mesh with a mesh generator */
199c4762a1bSJed Brown   PetscBool  testOrientIF;                 /* Test for different original interface orientations */
200c4762a1bSJed Brown   PetscBool  testHeavy;                    /* Run the heavy PointSF test */
201c4762a1bSJed Brown   PetscBool  customView;                   /* Show results of DMPlexIsInterpolated() etc. */
202c4762a1bSJed Brown   PetscInt   ornt[2];                      /* Orientation of interface on rank 0 and rank 1 */
203c4762a1bSJed Brown   PetscInt   faces[3];                     /* Number of faces per dimension for generator */
204d3e1f4ccSVaclav Hapla   PetscScalar coords[128];
205c4762a1bSJed Brown   PetscReal  coordsTol;
206c4762a1bSJed Brown   PetscInt   ncoords;
207c4762a1bSJed Brown   PetscInt   pointsToExpand[128];
208c4762a1bSJed Brown   PetscInt   nPointsToExpand;
209c4762a1bSJed Brown   PetscBool  testExpandPointsEmpty;
210c4762a1bSJed Brown   char       filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
211c4762a1bSJed Brown } AppCtx;
212c4762a1bSJed Brown 
213c4762a1bSJed Brown struct _n_PortableBoundary {
214c4762a1bSJed Brown   Vec coordinates;
215c4762a1bSJed Brown   PetscInt depth;
216c4762a1bSJed Brown   PetscSection *sections;
217c4762a1bSJed Brown };
218c4762a1bSJed Brown typedef struct _n_PortableBoundary * PortableBoundary;
219c4762a1bSJed Brown 
220956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
221c4762a1bSJed Brown static PetscLogStage  stage[3];
222956f8c0dSBarry Smith #endif
223c4762a1bSJed Brown 
224c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM, PortableBoundary);
225c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM,PetscBool);
226c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM, PortableBoundary *);
227c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM, IS, PetscSection, IS *);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown static PetscErrorCode PortableBoundaryDestroy(PortableBoundary *bnd)
230c4762a1bSJed Brown {
231c4762a1bSJed Brown   PetscInt       d;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   PetscFunctionBegin;
234c4762a1bSJed Brown   if (!*bnd) PetscFunctionReturn(0);
2359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*bnd)->coordinates));
236c4762a1bSJed Brown   for (d=0; d < (*bnd)->depth; d++) {
2379566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&(*bnd)->sections[d]));
238c4762a1bSJed Brown   }
2399566063dSJacob Faibussowitsch   PetscCall(PetscFree((*bnd)->sections));
2409566063dSJacob Faibussowitsch   PetscCall(PetscFree(*bnd));
241c4762a1bSJed Brown   PetscFunctionReturn(0);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
245c4762a1bSJed Brown {
246c4762a1bSJed Brown   const char    *interpTypes[4]  = {"none", "create", "after_create", "after_distribute"};
247c4762a1bSJed Brown   PetscInt       interp=NONE, dim;
248c4762a1bSJed Brown   PetscBool      flg1, flg2;
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   PetscFunctionBegin;
251c4762a1bSJed Brown   options->debug        = 0;
252c4762a1bSJed Brown   options->testNum      = 0;
253c4762a1bSJed Brown   options->dim          = 2;
254c4762a1bSJed Brown   options->cellSimplex  = PETSC_TRUE;
255c4762a1bSJed Brown   options->distribute   = PETSC_FALSE;
256c4762a1bSJed Brown   options->interpolate  = NONE;
257c4762a1bSJed Brown   options->useGenerator = PETSC_FALSE;
258c4762a1bSJed Brown   options->testOrientIF = PETSC_FALSE;
259c4762a1bSJed Brown   options->testHeavy    = PETSC_TRUE;
260c4762a1bSJed Brown   options->customView   = PETSC_FALSE;
261c4762a1bSJed Brown   options->testExpandPointsEmpty = PETSC_FALSE;
262c4762a1bSJed Brown   options->ornt[0]      = 0;
263c4762a1bSJed Brown   options->ornt[1]      = 0;
264c4762a1bSJed Brown   options->faces[0]     = 2;
265c4762a1bSJed Brown   options->faces[1]     = 2;
266c4762a1bSJed Brown   options->faces[2]     = 2;
267c4762a1bSJed Brown   options->filename[0]  = '\0';
268c4762a1bSJed Brown   options->coordsTol    = PETSC_DEFAULT;
269c4762a1bSJed Brown 
270d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
2719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex18.c", options->debug, &options->debug, NULL,0));
2729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex18.c", options->testNum, &options->testNum, NULL,0));
2739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_simplex", "Generate simplices if true, otherwise hexes", "ex18.c", options->cellSimplex, &options->cellSimplex, NULL));
2749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-distribute", "Distribute the mesh", "ex18.c", options->distribute, &options->distribute, NULL));
2759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEList("-interpolate", "Type of mesh interpolation (none, create, after_create, after_distribute)", "ex18.c", interpTypes, 4, interpTypes[options->interpolate], &interp, NULL));
276c4762a1bSJed Brown   options->interpolate = (InterpType) interp;
277*1dca8a05SBarry Smith   PetscCheck(options->distribute || options->interpolate != AFTER_DISTRIBUTE,comm, PETSC_ERR_SUP, "-interpolate after_distribute  needs  -distribute 1");
2789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex18.c", options->useGenerator, &options->useGenerator, NULL));
279c4762a1bSJed Brown   options->ncoords = 128;
2809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsScalarArray("-view_vertices_from_coords", "Print DAG points corresponding to vertices with given coordinates", "ex18.c", options->coords, &options->ncoords, NULL));
2819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-view_vertices_from_coords_tol", "Tolerance for -view_vertices_from_coords", "ex18.c", options->coordsTol, &options->coordsTol, NULL));
282c4762a1bSJed Brown   options->nPointsToExpand = 128;
2839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-test_expand_points", "Expand given array of DAG point using DMPlexGetConeRecursive() and print results", "ex18.c", options->pointsToExpand, &options->nPointsToExpand, NULL));
284c4762a1bSJed Brown   if (options->nPointsToExpand) {
2859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-test_expand_points_empty", "For -test_expand_points, rank 0 will have empty input array", "ex18.c", options->testExpandPointsEmpty, &options->testExpandPointsEmpty, NULL));
286c4762a1bSJed Brown   }
2879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_heavy", "Run the heavy PointSF test", "ex18.c", options->testHeavy, &options->testHeavy, NULL));
2889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-custom_view", "Custom DMPlex view", "ex18.c", options->customView, &options->customView, NULL));
2899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex18.c", options->dim, &options->dim, &flg1,1,3));
290c4762a1bSJed Brown   dim = 3;
2919566063dSJacob Faibussowitsch   PetscCall(PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex18.c", options->faces, &dim, &flg2));
292c4762a1bSJed Brown   if (flg2) {
293*1dca8a05SBarry Smith     PetscCheck(!flg1 || dim == options->dim,comm, PETSC_ERR_ARG_OUTOFRANGE, "specified -dim %" PetscInt_FMT " is not equal to length %" PetscInt_FMT " of -faces (note that -dim can be omitted)", options->dim, dim);
294c4762a1bSJed Brown     options->dim = dim;
295c4762a1bSJed Brown   }
2969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-filename", "The mesh file", "ex18.c", options->filename, options->filename, sizeof(options->filename), NULL));
2979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-rotate_interface_0", "Rotation (relative orientation) of interface on rank 0; implies -interpolate create -distribute 0", "ex18.c", options->ornt[0], &options->ornt[0], &options->testOrientIF,0));
2989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-rotate_interface_1", "Rotation (relative orientation) of interface on rank 1; implies -interpolate create -distribute 0", "ex18.c", options->ornt[1], &options->ornt[1], &flg2,0));
29908401ef6SPierre Jolivet   PetscCheck(flg2 == options->testOrientIF,comm, PETSC_ERR_ARG_OUTOFRANGE, "neither or both -rotate_interface_0 -rotate_interface_1 must be set");
300c4762a1bSJed Brown   if (options->testOrientIF) {
301c4762a1bSJed Brown     PetscInt i;
302c4762a1bSJed Brown     for (i=0; i<2; i++) {
303c4762a1bSJed Brown       if (options->ornt[i] >= 10) options->ornt[i] = -(options->ornt[i]-10);  /* 11 12 13 become -1 -2 -3 */
304c4762a1bSJed Brown     }
305c4762a1bSJed Brown     options->filename[0]  = 0;
306c4762a1bSJed Brown     options->useGenerator = PETSC_FALSE;
307c4762a1bSJed Brown     options->dim          = 3;
308c4762a1bSJed Brown     options->cellSimplex  = PETSC_TRUE;
309c4762a1bSJed Brown     options->interpolate  = CREATE;
310c4762a1bSJed Brown     options->distribute   = PETSC_FALSE;
311c4762a1bSJed Brown   }
312d0609cedSBarry Smith   PetscOptionsEnd();
313c4762a1bSJed Brown   PetscFunctionReturn(0);
314c4762a1bSJed Brown }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown static PetscErrorCode CreateMesh_1D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
317c4762a1bSJed Brown {
318c4762a1bSJed Brown   PetscInt       testNum = user->testNum;
319c4762a1bSJed Brown   PetscMPIInt    rank,size;
320c4762a1bSJed Brown   PetscInt       numCorners=2,i;
321c4762a1bSJed Brown   PetscInt       numCells,numVertices,network;
322a4a685f2SJacob Faibussowitsch   PetscInt       *cells;
323c4762a1bSJed Brown   PetscReal      *coords;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown   PetscFunctionBegin;
3269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
32863a3b9bcSJacob Faibussowitsch   PetscCheck(size <= 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for <=2 processes",testNum);
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   numCells = 3;
3319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-ncells", &numCells, NULL));
33263a3b9bcSJacob Faibussowitsch   PetscCheck(numCells >= 3,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test ncells %" PetscInt_FMT " must >=3",numCells);
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   if (size == 1) {
335c4762a1bSJed Brown     numVertices = numCells + 1;
3369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(2*numCells,&cells,2*numVertices,&coords));
337c4762a1bSJed Brown     for (i=0; i<numCells; i++) {
338c4762a1bSJed Brown       cells[2*i] = i; cells[2*i+1] = i + 1;
33971bbd31fSVaclav Hapla       coords[2*i] = i; coords[2*i+1] = i + 1;
340c4762a1bSJed Brown     }
341c4762a1bSJed Brown 
3429566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromCellListPetsc(comm, user->dim, numCells, numVertices, numCorners, PETSC_FALSE, cells, user->dim, coords, dm));
3439566063dSJacob Faibussowitsch     PetscCall(PetscFree2(cells,coords));
344c4762a1bSJed Brown     PetscFunctionReturn(0);
345c4762a1bSJed Brown   }
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   network = 0;
3489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-network_case", &network, NULL));
349c4762a1bSJed Brown   if (network == 0) {
350c4762a1bSJed Brown     switch (rank) {
351c4762a1bSJed Brown     case 0:
352c4762a1bSJed Brown     {
353c4762a1bSJed Brown       numCells    = 2;
354c4762a1bSJed Brown       numVertices = numCells;
3559566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(2*numCells,&cells,2*numCells,&coords));
356c4762a1bSJed Brown       cells[0] = 0; cells[1] = 1;
357c4762a1bSJed Brown       cells[2] = 1; cells[3] = 2;
358c4762a1bSJed Brown       coords[0] = 0.; coords[1] = 1.;
359c4762a1bSJed Brown       coords[2] = 1.; coords[3] = 2.;
360c4762a1bSJed Brown     }
361c4762a1bSJed Brown     break;
362c4762a1bSJed Brown     case 1:
363c4762a1bSJed Brown     {
364c4762a1bSJed Brown       numCells    -= 2;
365c4762a1bSJed Brown       numVertices = numCells + 1;
3669566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(2*numCells,&cells,2*numCells,&coords));
367c4762a1bSJed Brown       for (i=0; i<numCells; i++) {
368c4762a1bSJed Brown         cells[2*i] = 2+i; cells[2*i+1] = 2 + i + 1;
369c4762a1bSJed Brown         coords[2*i] = 2+i; coords[2*i+1] = 2 + i + 1;
370c4762a1bSJed Brown       }
371c4762a1bSJed Brown     }
372c4762a1bSJed Brown     break;
37398921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
374c4762a1bSJed Brown     }
375c4762a1bSJed Brown   } else { /* network_case = 1 */
376c4762a1bSJed Brown     /* ----------------------- */
377c4762a1bSJed Brown     switch (rank) {
378c4762a1bSJed Brown     case 0:
379c4762a1bSJed Brown     {
380c4762a1bSJed Brown       numCells    = 2;
381c4762a1bSJed Brown       numVertices = 3;
3829566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(2*numCells,&cells,2*numCells,&coords));
383c4762a1bSJed Brown       cells[0] = 0; cells[1] = 3;
384c4762a1bSJed Brown       cells[2] = 3; cells[3] = 1;
385c4762a1bSJed Brown     }
386c4762a1bSJed Brown     break;
387c4762a1bSJed Brown     case 1:
388c4762a1bSJed Brown     {
389c4762a1bSJed Brown       numCells    = 1;
390c4762a1bSJed Brown       numVertices = 1;
3919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(2*numCells,&cells,2*numCells,&coords));
392c4762a1bSJed Brown       cells[0] = 3; cells[1] = 2;
393c4762a1bSJed Brown     }
394c4762a1bSJed Brown     break;
39598921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
396c4762a1bSJed Brown     }
397c4762a1bSJed Brown   }
3989566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, PETSC_FALSE, cells, user->dim, coords, NULL, NULL, dm));
3999566063dSJacob Faibussowitsch   PetscCall(PetscFree2(cells,coords));
400c4762a1bSJed Brown   PetscFunctionReturn(0);
401c4762a1bSJed Brown }
402c4762a1bSJed Brown 
403c4762a1bSJed Brown static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
404c4762a1bSJed Brown {
405c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
406c4762a1bSJed Brown   PetscMPIInt    rank, size;
407c4762a1bSJed Brown 
408c4762a1bSJed Brown   PetscFunctionBegin;
4099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
411c4762a1bSJed Brown   switch (testNum) {
412c4762a1bSJed Brown   case 0:
41363a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum);
414c4762a1bSJed Brown     switch (rank) {
415c4762a1bSJed Brown       case 0:
416c4762a1bSJed Brown       {
417c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
418a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {0, 1, 2};
419c4762a1bSJed Brown         PetscReal      coords[4] = {-0.5, 0.5, 0.0, 0.0};
420c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
421c4762a1bSJed Brown 
4229566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4239566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
424c4762a1bSJed Brown       }
425c4762a1bSJed Brown       break;
426c4762a1bSJed Brown       case 1:
427c4762a1bSJed Brown       {
428c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
429a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {1, 3, 2};
430c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 1.0, 0.5, 0.5};
431c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
432c4762a1bSJed Brown 
4339566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4349566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
435c4762a1bSJed Brown       }
436c4762a1bSJed Brown       break;
43798921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
438c4762a1bSJed Brown     }
439c4762a1bSJed Brown     break;
440c4762a1bSJed Brown   case 1:
44163a3b9bcSJacob Faibussowitsch     PetscCheck(size == 3,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum);
442c4762a1bSJed Brown     switch (rank) {
443c4762a1bSJed Brown       case 0:
444c4762a1bSJed Brown       {
445c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
446a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {0, 1, 2};
447c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 1.0, 0.0, 0.0};
448c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
449c4762a1bSJed Brown 
4509566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4519566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
452c4762a1bSJed Brown       }
453c4762a1bSJed Brown       break;
454c4762a1bSJed Brown       case 1:
455c4762a1bSJed Brown       {
456c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
457a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {0, 2, 3};
458c4762a1bSJed Brown         PetscReal      coords[4] = {0.5, 0.5, 1.0, 1.0};
459c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
460c4762a1bSJed Brown 
4619566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4629566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
463c4762a1bSJed Brown       }
464c4762a1bSJed Brown       break;
465c4762a1bSJed Brown       case 2:
466c4762a1bSJed Brown       {
467c4762a1bSJed Brown         const PetscInt numCells  = 2, numVertices = 1, numCorners = 3;
468a4a685f2SJacob Faibussowitsch         const PetscInt cells[6]  = {2, 4, 3, 2, 1, 4};
469c4762a1bSJed Brown         PetscReal      coords[2] = {1.0, 0.0};
470c4762a1bSJed Brown         PetscInt       markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
471c4762a1bSJed Brown 
4729566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4739566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
474c4762a1bSJed Brown       }
475c4762a1bSJed Brown       break;
47698921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
477c4762a1bSJed Brown     }
478c4762a1bSJed Brown     break;
479c4762a1bSJed Brown   case 2:
48063a3b9bcSJacob Faibussowitsch     PetscCheck(size == 3,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 3 processes", testNum);
481c4762a1bSJed Brown     switch (rank) {
482c4762a1bSJed Brown       case 0:
483c4762a1bSJed Brown       {
484c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
485a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {1, 2, 0};
486c4762a1bSJed Brown         PetscReal      coords[4] = {0.5, 0.5, 0.0, 1.0};
487c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
488c4762a1bSJed Brown 
4899566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
4909566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
491c4762a1bSJed Brown       }
492c4762a1bSJed Brown       break;
493c4762a1bSJed Brown       case 1:
494c4762a1bSJed Brown       {
495c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
496a4a685f2SJacob Faibussowitsch         const PetscInt cells[3]  = {1, 0, 3};
497c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 0.0, 1.0, 1.0};
498c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
499c4762a1bSJed Brown 
5009566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
5019566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
502c4762a1bSJed Brown       }
503c4762a1bSJed Brown       break;
504c4762a1bSJed Brown       case 2:
505c4762a1bSJed Brown       {
506c4762a1bSJed Brown         const PetscInt numCells  = 2, numVertices = 1, numCorners = 3;
507a4a685f2SJacob Faibussowitsch         const PetscInt cells[6]  = {0, 4, 3, 0, 2, 4};
508c4762a1bSJed Brown         PetscReal      coords[2] = {1.0, 0.0};
509c4762a1bSJed Brown         PetscInt       markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
510c4762a1bSJed Brown 
5119566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
5129566063dSJacob Faibussowitsch         for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
513c4762a1bSJed Brown       }
514c4762a1bSJed Brown       break;
51598921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
516c4762a1bSJed Brown     }
517c4762a1bSJed Brown     break;
51863a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
519c4762a1bSJed Brown   }
520c4762a1bSJed Brown   PetscFunctionReturn(0);
521c4762a1bSJed Brown }
522c4762a1bSJed Brown 
523c4762a1bSJed Brown static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
524c4762a1bSJed Brown {
525c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
526c4762a1bSJed Brown   PetscMPIInt    rank, size;
527c4762a1bSJed Brown 
528c4762a1bSJed Brown   PetscFunctionBegin;
5299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
531c4762a1bSJed Brown   switch (testNum) {
532c4762a1bSJed Brown   case 0:
53363a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum);
534c4762a1bSJed Brown     switch (rank) {
535c4762a1bSJed Brown       case 0:
536c4762a1bSJed Brown       {
537c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 4;
538a4a685f2SJacob Faibussowitsch         const PetscInt cells[4]  = {0, 2, 1, 3};
539c4762a1bSJed Brown         PetscReal      coords[6] = {0.0, 0.0, -0.5,  0.0, -0.5, 0.0};
540c4762a1bSJed Brown         PetscInt       markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1};
541c4762a1bSJed Brown 
5429566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
5439566063dSJacob Faibussowitsch         for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
544c4762a1bSJed Brown       }
545c4762a1bSJed Brown       break;
546c4762a1bSJed Brown       case 1:
547c4762a1bSJed Brown       {
548c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
549a4a685f2SJacob Faibussowitsch         const PetscInt cells[4]  = {1, 2, 4, 3};
550c4762a1bSJed Brown         PetscReal      coords[9] = {1.0, 0.0, 0.0,  0.0, 0.5, 0.0,  0.0, 0.0, 0.5};
551c4762a1bSJed Brown         PetscInt       markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1};
552c4762a1bSJed Brown 
5539566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
5549566063dSJacob Faibussowitsch         for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
555c4762a1bSJed Brown       }
556c4762a1bSJed Brown       break;
55798921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
558c4762a1bSJed Brown     }
559c4762a1bSJed Brown     break;
56063a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
561c4762a1bSJed Brown   }
562c4762a1bSJed Brown   if (user->testOrientIF) {
563c4762a1bSJed Brown     PetscInt ifp[] = {8, 6};
564c4762a1bSJed Brown 
5659566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) *dm, "Mesh before orientation"));
5669566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view"));
567c4762a1bSJed Brown     /* rotate interface face ifp[rank] by given orientation ornt[rank] */
5689566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientPoint(*dm, ifp[rank], user->ornt[rank]));
5699566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view"));
5709566063dSJacob Faibussowitsch     PetscCall(DMPlexCheckFaces(*dm, 0));
5719566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientInterface_Internal(*dm));
5729566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Orientation test PASSED\n"));
573c4762a1bSJed Brown   }
574c4762a1bSJed Brown   PetscFunctionReturn(0);
575c4762a1bSJed Brown }
576c4762a1bSJed Brown 
577c4762a1bSJed Brown static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
578c4762a1bSJed Brown {
579c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
580c4762a1bSJed Brown   PetscMPIInt    rank, size;
581c4762a1bSJed Brown 
582c4762a1bSJed Brown   PetscFunctionBegin;
5839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
585c4762a1bSJed Brown   switch (testNum) {
586c4762a1bSJed Brown   case 0:
58763a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum);
588c4762a1bSJed Brown     switch (rank) {
589c4762a1bSJed Brown       case 0:
590c4762a1bSJed Brown       {
591c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
592a4a685f2SJacob Faibussowitsch         const PetscInt cells[4]  = {0, 1, 2, 3};
593c4762a1bSJed Brown         PetscReal      coords[6] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0};
594c4762a1bSJed Brown         PetscInt       markerPoints[4*2] = {1, 1, 2, 1, 3, 1, 4, 1};
595c4762a1bSJed Brown 
5969566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
5979566063dSJacob Faibussowitsch         for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
598c4762a1bSJed Brown       }
599c4762a1bSJed Brown       break;
600c4762a1bSJed Brown       case 1:
601c4762a1bSJed Brown       {
602c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
603a4a685f2SJacob Faibussowitsch         const PetscInt cells[4]  = {1, 4, 5, 2};
604c4762a1bSJed Brown         PetscReal      coords[6] = {-0.5, 1.0, 0.5, 0.0, 0.5, 1.0};
605c4762a1bSJed Brown         PetscInt       markerPoints[4*2] = {1, 1, 2, 1, 3, 1, 4, 1};
606c4762a1bSJed Brown 
6079566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
6089566063dSJacob Faibussowitsch         for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
609c4762a1bSJed Brown       }
610c4762a1bSJed Brown       break;
61198921bdaSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
612c4762a1bSJed Brown     }
613c4762a1bSJed Brown     break;
61463a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
615c4762a1bSJed Brown   }
616c4762a1bSJed Brown   PetscFunctionReturn(0);
617c4762a1bSJed Brown }
618c4762a1bSJed Brown 
619c4762a1bSJed Brown static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
620c4762a1bSJed Brown {
621c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
622c4762a1bSJed Brown   PetscMPIInt    rank, size;
623c4762a1bSJed Brown 
624c4762a1bSJed Brown   PetscFunctionBegin;
6259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
6269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
627c4762a1bSJed Brown   switch (testNum) {
628c4762a1bSJed Brown   case 0:
62963a3b9bcSJacob Faibussowitsch     PetscCheck(size == 2,comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %" PetscInt_FMT " only for 2 processes", testNum);
630c4762a1bSJed Brown     switch (rank) {
631c4762a1bSJed Brown     case 0:
632c4762a1bSJed Brown     {
633c4762a1bSJed Brown       const PetscInt numCells  = 1, numVertices = 6, numCorners = 8;
634a4a685f2SJacob Faibussowitsch       const PetscInt cells[8]  = {0, 3, 2, 1, 4, 5, 6, 7};
635c4762a1bSJed Brown       PetscReal      coords[6*3] = {-0.5,0.0,0.0, 0.0,0.0,0.0, 0.0,1.0,0.0, -0.5,1.0,0.0, -0.5,0.0,1.0, 0.0,0.0,1.0};
636c4762a1bSJed Brown       PetscInt       markerPoints[8*2] = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
637c4762a1bSJed Brown 
6389566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
6399566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
640c4762a1bSJed Brown     }
641c4762a1bSJed Brown     break;
642c4762a1bSJed Brown     case 1:
643c4762a1bSJed Brown     {
644c4762a1bSJed Brown       const PetscInt numCells  = 1, numVertices = 6, numCorners = 8;
645a4a685f2SJacob Faibussowitsch       const PetscInt cells[8]  = {1, 2, 9, 8, 5, 10, 11, 6};
646c4762a1bSJed Brown       PetscReal      coords[6*3] = {0.0,1.0,1.0, -0.5,1.0,1.0, 0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0};
647c4762a1bSJed Brown       PetscInt       markerPoints[8*2] = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
648c4762a1bSJed Brown 
6499566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromCellListParallelPetsc(comm, user->dim, numCells, numVertices, PETSC_DECIDE, numCorners, interpolate, cells, user->dim, coords, NULL, NULL, dm));
6509566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
651c4762a1bSJed Brown     }
652c4762a1bSJed Brown     break;
65398921bdaSJacob Faibussowitsch     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
654c4762a1bSJed Brown     }
655c4762a1bSJed Brown   break;
65663a3b9bcSJacob Faibussowitsch   default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
657c4762a1bSJed Brown   }
658c4762a1bSJed Brown   PetscFunctionReturn(0);
659c4762a1bSJed Brown }
660c4762a1bSJed Brown 
661c4762a1bSJed Brown static PetscErrorCode CustomView(DM dm, PetscViewer v)
662c4762a1bSJed Brown {
663c4762a1bSJed Brown   DMPlexInterpolatedFlag interpolated;
664c4762a1bSJed Brown   PetscBool              distributed;
665c4762a1bSJed Brown 
666c4762a1bSJed Brown   PetscFunctionBegin;
6679566063dSJacob Faibussowitsch   PetscCall(DMPlexIsDistributed(dm, &distributed));
6689566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(dm, &interpolated));
6699566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsDistributed: %s\n", PetscBools[distributed]));
6709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(v, "DMPlexIsInterpolatedCollective: %s\n", DMPlexInterpolatedFlags[interpolated]));
671c4762a1bSJed Brown   PetscFunctionReturn(0);
672c4762a1bSJed Brown }
673c4762a1bSJed Brown 
674c4762a1bSJed Brown static PetscErrorCode CreateMeshFromFile(MPI_Comm comm, AppCtx *user, DM *dm, DM *serialDM)
675c4762a1bSJed Brown {
676c4762a1bSJed Brown   const char    *filename       = user->filename;
677c4762a1bSJed Brown   PetscBool      testHeavy      = user->testHeavy;
678c4762a1bSJed Brown   PetscBool      interpCreate   = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE;
679c4762a1bSJed Brown   PetscBool      distributed    = PETSC_FALSE;
680c4762a1bSJed Brown 
681c4762a1bSJed Brown   PetscFunctionBegin;
682c4762a1bSJed Brown   *serialDM = NULL;
6839566063dSJacob Faibussowitsch   if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_FALSE));
6849566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage[0]));
6859566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, dm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
6869566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
6879566063dSJacob Faibussowitsch   if (testHeavy && interpCreate) PetscCall(DMPlexSetOrientInterface_Private(NULL, PETSC_TRUE));
6889566063dSJacob Faibussowitsch   PetscCall(DMPlexIsDistributed(*dm, &distributed));
6899566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "DMPlexCreateFromFile produced %s mesh.\n", distributed ? "distributed" : "serial"));
690c4762a1bSJed Brown   if (testHeavy && distributed) {
6919566063dSJacob Faibussowitsch     PetscCall(PetscOptionsSetValue(NULL, "-dm_plex_hdf5_force_sequential", NULL));
6929566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromFile(comm, filename, "ex18_plex", interpCreate, serialDM));
6939566063dSJacob Faibussowitsch     PetscCall(DMPlexIsDistributed(*serialDM, &distributed));
69428b400f6SJacob Faibussowitsch     PetscCheck(!distributed,comm, PETSC_ERR_PLIB, "unable to create a serial DM from file");
695c4762a1bSJed Brown   }
6969566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(*dm, &user->dim));
697c4762a1bSJed Brown   PetscFunctionReturn(0);
698c4762a1bSJed Brown }
699c4762a1bSJed Brown 
700c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
701c4762a1bSJed Brown {
702c4762a1bSJed Brown   PetscPartitioner part;
703c4762a1bSJed Brown   PortableBoundary boundary     = NULL;
704c4762a1bSJed Brown   DM             serialDM       = NULL;
705c4762a1bSJed Brown   PetscBool      cellSimplex    = user->cellSimplex;
706c4762a1bSJed Brown   PetscBool      useGenerator   = user->useGenerator;
707c4762a1bSJed Brown   PetscBool      interpCreate   = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE;
708c4762a1bSJed Brown   PetscBool      interpSerial   = user->interpolate == AFTER_CREATE ? PETSC_TRUE : PETSC_FALSE;
709c4762a1bSJed Brown   PetscBool      interpParallel = user->interpolate == AFTER_DISTRIBUTE ? PETSC_TRUE : PETSC_FALSE;
710c4762a1bSJed Brown   PetscBool      testHeavy      = user->testHeavy;
711c4762a1bSJed Brown   PetscMPIInt    rank;
712c4762a1bSJed Brown 
713c4762a1bSJed Brown   PetscFunctionBegin;
7149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
715c4762a1bSJed Brown   if (user->filename[0]) {
7169566063dSJacob Faibussowitsch     PetscCall(CreateMeshFromFile(comm, user, dm, &serialDM));
717c4762a1bSJed Brown   } else if (useGenerator) {
7189566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[0]));
7199566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateBoxMesh(comm, user->dim, cellSimplex, user->faces, NULL, NULL, NULL, interpCreate, dm));
7209566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
721c4762a1bSJed Brown   } else {
7229566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[0]));
723c4762a1bSJed Brown     switch (user->dim) {
724c4762a1bSJed Brown     case 1:
7259566063dSJacob Faibussowitsch       PetscCall(CreateMesh_1D(comm, interpCreate, user, dm));
726c4762a1bSJed Brown       break;
727c4762a1bSJed Brown     case 2:
728c4762a1bSJed Brown       if (cellSimplex) {
7299566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_2D(comm, interpCreate, user, dm));
730c4762a1bSJed Brown       } else {
7319566063dSJacob Faibussowitsch         PetscCall(CreateQuad_2D(comm, interpCreate, user, dm));
732c4762a1bSJed Brown       }
733c4762a1bSJed Brown       break;
734c4762a1bSJed Brown     case 3:
735c4762a1bSJed Brown       if (cellSimplex) {
7369566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_3D(comm, interpCreate, user, dm));
737c4762a1bSJed Brown       } else {
7389566063dSJacob Faibussowitsch         PetscCall(CreateHex_3D(comm, interpCreate, user, dm));
739c4762a1bSJed Brown       }
740c4762a1bSJed Brown       break;
741c4762a1bSJed Brown     default:
74263a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, user->dim);
743c4762a1bSJed Brown     }
7449566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
745c4762a1bSJed Brown   }
74663a3b9bcSJacob Faibussowitsch   PetscCheck(user->ncoords % user->dim == 0,comm, PETSC_ERR_ARG_OUTOFRANGE, "length of coordinates array %" PetscInt_FMT " must be divisable by spatial dimension %" PetscInt_FMT, user->ncoords, user->dim);
7479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "Original Mesh"));
7489566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-orig_dm_view"));
749c4762a1bSJed Brown 
750c4762a1bSJed Brown   if (interpSerial) {
751c4762a1bSJed Brown     DM idm;
752c4762a1bSJed Brown 
7539566063dSJacob Faibussowitsch     if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE));
7549566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[2]));
7559566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
7569566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
7579566063dSJacob Faibussowitsch     if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE));
7589566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
759c4762a1bSJed Brown     *dm = idm;
7609566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) *dm, "Interpolated Mesh"));
7619566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view"));
762c4762a1bSJed Brown   }
763c4762a1bSJed Brown 
764c4762a1bSJed Brown   /* Set partitioner options */
7659566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(*dm, &part));
766c4762a1bSJed Brown   if (part) {
7679566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE));
7689566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
769c4762a1bSJed Brown   }
770c4762a1bSJed Brown 
7719566063dSJacob Faibussowitsch   if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm)));
772c4762a1bSJed Brown   if (testHeavy) {
773c4762a1bSJed Brown     PetscBool distributed;
774c4762a1bSJed Brown 
7759566063dSJacob Faibussowitsch     PetscCall(DMPlexIsDistributed(*dm, &distributed));
776c4762a1bSJed Brown     if (!serialDM && !distributed) {
777c4762a1bSJed Brown       serialDM = *dm;
7789566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)*dm));
779c4762a1bSJed Brown     }
780c4762a1bSJed Brown     if (serialDM) {
7819566063dSJacob Faibussowitsch       PetscCall(DMPlexGetExpandedBoundary_Private(serialDM, &boundary));
782c4762a1bSJed Brown     }
783c4762a1bSJed Brown     if (boundary) {
784c4762a1bSJed Brown       /* check DM which has been created in parallel and already interpolated */
7859566063dSJacob Faibussowitsch       PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary));
786c4762a1bSJed Brown     }
787c4762a1bSJed Brown     /* Orient interface because it could be deliberately skipped above. It is idempotent. */
7889566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientInterface_Internal(*dm));
789c4762a1bSJed Brown   }
790c4762a1bSJed Brown   if (user->distribute) {
791c4762a1bSJed Brown     DM               pdm = NULL;
792c4762a1bSJed Brown 
793c4762a1bSJed Brown     /* Redistribute mesh over processes using that partitioner */
7949566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage[1]));
7959566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
7969566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
797c4762a1bSJed Brown     if (pdm) {
7989566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
799c4762a1bSJed Brown       *dm  = pdm;
8009566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) *dm, "Redistributed Mesh"));
8019566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(*dm, NULL, "-dist_dm_view"));
802c4762a1bSJed Brown     }
803c4762a1bSJed Brown 
804c4762a1bSJed Brown     if (interpParallel) {
805c4762a1bSJed Brown       DM idm;
806c4762a1bSJed Brown 
8079566063dSJacob Faibussowitsch       if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE));
8089566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePush(stage[2]));
8099566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
8109566063dSJacob Faibussowitsch       PetscCall(PetscLogStagePop());
8119566063dSJacob Faibussowitsch       if (testHeavy) PetscCall(DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE));
8129566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
813c4762a1bSJed Brown       *dm = idm;
8149566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) *dm, "Interpolated Redistributed Mesh"));
8159566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(*dm, NULL, "-intp_dm_view"));
816c4762a1bSJed Brown     }
817c4762a1bSJed Brown   }
818c4762a1bSJed Brown   if (testHeavy) {
819c4762a1bSJed Brown     if (boundary) {
8209566063dSJacob Faibussowitsch       PetscCall(DMPlexCheckPointSFHeavy(*dm, boundary));
821c4762a1bSJed Brown     }
822c4762a1bSJed Brown     /* Orient interface because it could be deliberately skipped above. It is idempotent. */
8239566063dSJacob Faibussowitsch     PetscCall(DMPlexOrientInterface_Internal(*dm));
824c4762a1bSJed Brown   }
825c4762a1bSJed Brown 
8269566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "Parallel Mesh"));
8279566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
8289566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
8299566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
830c4762a1bSJed Brown 
8319566063dSJacob Faibussowitsch   if (user->customView) PetscCall(CustomView(*dm, PETSC_VIEWER_STDOUT_(comm)));
8329566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&serialDM));
8339566063dSJacob Faibussowitsch   PetscCall(PortableBoundaryDestroy(&boundary));
834c4762a1bSJed Brown   PetscFunctionReturn(0);
835c4762a1bSJed Brown }
836c4762a1bSJed Brown 
837d3e1f4ccSVaclav Hapla #define ps2d(number) ((double) PetscRealPart(number))
8389fbee547SJacob Faibussowitsch static inline PetscErrorCode coord2str(char buf[], size_t len, PetscInt dim, const PetscScalar coords[], PetscReal tol)
839c4762a1bSJed Brown {
840c4762a1bSJed Brown   PetscFunctionBegin;
84108401ef6SPierre Jolivet   PetscCheck(dim <= 3,PETSC_COMM_SELF, PETSC_ERR_SUP, "dim must be less than or equal 3");
842c4762a1bSJed Brown   if (tol >= 1e-3) {
843c4762a1bSJed Brown     switch (dim) {
8449566063dSJacob Faibussowitsch       case 1: PetscCall(PetscSNPrintf(buf,len,"(%12.3f)",ps2d(coords[0])));
8459566063dSJacob Faibussowitsch       case 2: PetscCall(PetscSNPrintf(buf,len,"(%12.3f, %12.3f)",ps2d(coords[0]),ps2d(coords[1])));
8469566063dSJacob Faibussowitsch       case 3: PetscCall(PetscSNPrintf(buf,len,"(%12.3f, %12.3f, %12.3f)",ps2d(coords[0]),ps2d(coords[1]),ps2d(coords[2])));
847c4762a1bSJed Brown     }
848c4762a1bSJed Brown   } else {
849c4762a1bSJed Brown     switch (dim) {
8509566063dSJacob Faibussowitsch       case 1: PetscCall(PetscSNPrintf(buf,len,"(%12.6f)",ps2d(coords[0])));
8519566063dSJacob Faibussowitsch       case 2: PetscCall(PetscSNPrintf(buf,len,"(%12.6f, %12.6f)",ps2d(coords[0]),ps2d(coords[1])));
8529566063dSJacob Faibussowitsch       case 3: PetscCall(PetscSNPrintf(buf,len,"(%12.6f, %12.6f, %12.6f)",ps2d(coords[0]),ps2d(coords[1]),ps2d(coords[2])));
853c4762a1bSJed Brown     }
854c4762a1bSJed Brown   }
855c4762a1bSJed Brown   PetscFunctionReturn(0);
856c4762a1bSJed Brown }
857c4762a1bSJed Brown 
858d3e1f4ccSVaclav Hapla static PetscErrorCode ViewVerticesFromCoords(DM dm, Vec coordsVec, PetscReal tol, PetscViewer viewer)
859c4762a1bSJed Brown {
860d3e1f4ccSVaclav Hapla   PetscInt       dim, i, npoints;
861d3e1f4ccSVaclav Hapla   IS             pointsIS;
862d3e1f4ccSVaclav Hapla   const PetscInt *points;
863d3e1f4ccSVaclav Hapla   const PetscScalar *coords;
864c4762a1bSJed Brown   char           coordstr[128];
865c4762a1bSJed Brown   MPI_Comm       comm;
866c4762a1bSJed Brown   PetscMPIInt    rank;
867c4762a1bSJed Brown 
868c4762a1bSJed Brown   PetscFunctionBegin;
8699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
8709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8719566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
8729566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(viewer));
8739566063dSJacob Faibussowitsch   PetscCall(DMPlexFindVertices(dm, coordsVec, tol, &pointsIS));
8749566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointsIS, &points));
8759566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointsIS, &npoints));
8769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(coordsVec, &coords));
877c4762a1bSJed Brown   for (i=0; i < npoints; i++) {
8789566063dSJacob Faibussowitsch     PetscCall(coord2str(coordstr, sizeof(coordstr), dim, &coords[i*dim], tol));
8799566063dSJacob Faibussowitsch     if (rank == 0 && i) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "-----\n"));
88063a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %s --> points[%" PetscInt_FMT "] = %" PetscInt_FMT "\n", rank, coordstr, i, points[i]));
8819566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
882c4762a1bSJed Brown   }
8839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(viewer));
8849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(coordsVec, &coords));
8859566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointsIS, &points));
8869566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointsIS));
887c4762a1bSJed Brown   PetscFunctionReturn(0);
888c4762a1bSJed Brown }
889c4762a1bSJed Brown 
890c4762a1bSJed Brown static PetscErrorCode TestExpandPoints(DM dm, AppCtx *user)
891c4762a1bSJed Brown {
892c4762a1bSJed Brown   IS                is;
893c4762a1bSJed Brown   PetscSection      *sects;
894c4762a1bSJed Brown   IS                *iss;
895c4762a1bSJed Brown   PetscInt          d,depth;
896c4762a1bSJed Brown   PetscMPIInt       rank;
897c4762a1bSJed Brown   PetscViewer       viewer=PETSC_VIEWER_STDOUT_WORLD, sviewer;
898c4762a1bSJed Brown 
899c4762a1bSJed Brown   PetscFunctionBegin;
9009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
901dd400576SPatrick Sanan   if (user->testExpandPointsEmpty && rank == 0) {
9029566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &is));
903c4762a1bSJed Brown   } else {
9049566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_SELF, user->nPointsToExpand, user->pointsToExpand, PETSC_USE_POINTER, &is));
905c4762a1bSJed Brown   }
9069566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeRecursive(dm, is, &depth, &iss, &sects));
9079566063dSJacob Faibussowitsch   PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
9089566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] ==========================\n",rank));
909c4762a1bSJed Brown   for (d=depth-1; d>=0; d--) {
910c4762a1bSJed Brown     IS          checkIS;
911c4762a1bSJed Brown     PetscBool   flg;
912c4762a1bSJed Brown 
91363a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(sviewer, "depth %" PetscInt_FMT " ---------------\n",d));
9149566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sects[d], sviewer));
9159566063dSJacob Faibussowitsch     PetscCall(ISView(iss[d], sviewer));
916c4762a1bSJed Brown     /* check reverse operation */
917c4762a1bSJed Brown     if (d < depth-1) {
9189566063dSJacob Faibussowitsch       PetscCall(DMPlexExpandedConesToFaces_Private(dm, iss[d], sects[d], &checkIS));
9199566063dSJacob Faibussowitsch       PetscCall(ISEqualUnsorted(checkIS, iss[d+1], &flg));
92028b400f6SJacob Faibussowitsch       PetscCheck(flg,PetscObjectComm((PetscObject) checkIS), PETSC_ERR_PLIB, "DMPlexExpandedConesToFaces_Private produced wrong IS");
9219566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&checkIS));
922c4762a1bSJed Brown     }
923c4762a1bSJed Brown   }
9249566063dSJacob Faibussowitsch   PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer));
9259566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
9269566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreConeRecursive(dm, is, &depth, &iss, &sects));
9279566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is));
928c4762a1bSJed Brown   PetscFunctionReturn(0);
929c4762a1bSJed Brown }
930c4762a1bSJed Brown 
931c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM dm, IS is, PetscSection section, IS *newis)
932c4762a1bSJed Brown {
933c4762a1bSJed Brown   PetscInt          n,n1,ncone,numCoveredPoints,o,p,q,start,end;
934c4762a1bSJed Brown   const PetscInt    *coveredPoints;
935c4762a1bSJed Brown   const PetscInt    *arr, *cone;
936c4762a1bSJed Brown   PetscInt          *newarr;
937c4762a1bSJed Brown 
938c4762a1bSJed Brown   PetscFunctionBegin;
9399566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(is, &n));
9409566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &n1));
9419566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &start, &end));
94263a3b9bcSJacob Faibussowitsch   PetscCheck(n == n1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "IS size = %" PetscInt_FMT " != %" PetscInt_FMT " = section storage size", n, n1);
9439566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(is, &arr));
9449566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(end-start, &newarr));
945c4762a1bSJed Brown   for (q=start; q<end; q++) {
9469566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, q, &ncone));
9479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, q, &o));
948c4762a1bSJed Brown     cone = &arr[o];
949c4762a1bSJed Brown     if (ncone == 1) {
950c4762a1bSJed Brown       numCoveredPoints = 1;
951c4762a1bSJed Brown       p = cone[0];
952c4762a1bSJed Brown     } else {
953c4762a1bSJed Brown       PetscInt i;
954c4762a1bSJed Brown       p = PETSC_MAX_INT;
955c4762a1bSJed Brown       for (i=0; i<ncone; i++) if (cone[i] < 0) {p = -1; break;}
956c4762a1bSJed Brown       if (p >= 0) {
9579566063dSJacob Faibussowitsch         PetscCall(DMPlexGetJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints));
95863a3b9bcSJacob Faibussowitsch         PetscCheck(numCoveredPoints <= 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "more than one covered points for section point %" PetscInt_FMT,q);
959c4762a1bSJed Brown         if (numCoveredPoints) p = coveredPoints[0];
960c4762a1bSJed Brown         else                  p = -2;
9619566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints));
962c4762a1bSJed Brown       }
963c4762a1bSJed Brown     }
964c4762a1bSJed Brown     newarr[q-start] = p;
965c4762a1bSJed Brown   }
9669566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(is, &arr));
9679566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, end-start, newarr, PETSC_OWN_POINTER, newis));
968c4762a1bSJed Brown   PetscFunctionReturn(0);
969c4762a1bSJed Brown }
970c4762a1bSJed Brown 
971c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedVerticesToFaces_Private(DM dm, IS boundary_expanded_is, PetscInt depth, PetscSection sections[], IS *boundary_is)
972c4762a1bSJed Brown {
973c4762a1bSJed Brown   PetscInt          d;
974c4762a1bSJed Brown   IS                is,newis;
975c4762a1bSJed Brown 
976c4762a1bSJed Brown   PetscFunctionBegin;
977c4762a1bSJed Brown   is = boundary_expanded_is;
9789566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)is));
979c4762a1bSJed Brown   for (d = 0; d < depth-1; ++d) {
9809566063dSJacob Faibussowitsch     PetscCall(DMPlexExpandedConesToFaces_Private(dm, is, sections[d], &newis));
9819566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&is));
982c4762a1bSJed Brown     is = newis;
983c4762a1bSJed Brown   }
984c4762a1bSJed Brown   *boundary_is = is;
985c4762a1bSJed Brown   PetscFunctionReturn(0);
986c4762a1bSJed Brown }
987c4762a1bSJed Brown 
988c4762a1bSJed Brown #define CHKERRQI(incall,ierr) if (ierr) {incall = PETSC_FALSE; }
989c4762a1bSJed Brown 
990c4762a1bSJed Brown static PetscErrorCode DMLabelViewFromOptionsOnComm_Private(DMLabel label, const char optionname[], MPI_Comm comm)
991c4762a1bSJed Brown {
992c4762a1bSJed Brown   PetscViewer       viewer;
993c4762a1bSJed Brown   PetscBool         flg;
994c4762a1bSJed Brown   static PetscBool  incall = PETSC_FALSE;
995c4762a1bSJed Brown   PetscViewerFormat format;
996c4762a1bSJed Brown 
997c4762a1bSJed Brown   PetscFunctionBegin;
998c4762a1bSJed Brown   if (incall) PetscFunctionReturn(0);
999c4762a1bSJed Brown   incall = PETSC_TRUE;
10005f80ce2aSJacob Faibussowitsch   CHKERRQI(incall,PetscOptionsGetViewer(comm,((PetscObject)label)->options,((PetscObject)label)->prefix,optionname,&viewer,&format,&flg));
1001c4762a1bSJed Brown   if (flg) {
10025f80ce2aSJacob Faibussowitsch     CHKERRQI(incall,PetscViewerPushFormat(viewer,format));
10035f80ce2aSJacob Faibussowitsch     CHKERRQI(incall,DMLabelView(label, viewer));
10045f80ce2aSJacob Faibussowitsch     CHKERRQI(incall,PetscViewerPopFormat(viewer));
10055f80ce2aSJacob Faibussowitsch     CHKERRQI(incall,PetscViewerDestroy(&viewer));
1006c4762a1bSJed Brown   }
1007c4762a1bSJed Brown   incall = PETSC_FALSE;
1008c4762a1bSJed Brown   PetscFunctionReturn(0);
1009c4762a1bSJed Brown }
1010c4762a1bSJed Brown 
1011c4762a1bSJed Brown /* TODO: this is hotfixing DMLabelGetStratumIS() - it should be fixed systematically instead */
10129fbee547SJacob Faibussowitsch static inline PetscErrorCode DMLabelGetStratumISOnComm_Private(DMLabel label, PetscInt value, MPI_Comm comm, IS *is)
1013c4762a1bSJed Brown {
1014c4762a1bSJed Brown   IS                tmpis;
1015c4762a1bSJed Brown 
1016c4762a1bSJed Brown   PetscFunctionBegin;
10179566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(label, value, &tmpis));
10189566063dSJacob Faibussowitsch   if (!tmpis) PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &tmpis));
10199566063dSJacob Faibussowitsch   PetscCall(ISOnComm(tmpis, comm, PETSC_COPY_VALUES, is));
10209566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&tmpis));
1021c4762a1bSJed Brown   PetscFunctionReturn(0);
1022c4762a1bSJed Brown }
1023c4762a1bSJed Brown 
1024c4762a1bSJed Brown /* currently only for simple PetscSection without fields or constraints */
1025c4762a1bSJed Brown static PetscErrorCode PetscSectionReplicate_Private(MPI_Comm comm, PetscMPIInt rootrank, PetscSection sec0, PetscSection *secout)
1026c4762a1bSJed Brown {
1027c4762a1bSJed Brown   PetscSection      sec;
1028c4762a1bSJed Brown   PetscInt          chart[2], p;
1029c4762a1bSJed Brown   PetscInt          *dofarr;
1030c4762a1bSJed Brown   PetscMPIInt       rank;
1031c4762a1bSJed Brown 
1032c4762a1bSJed Brown   PetscFunctionBegin;
10339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1034c4762a1bSJed Brown   if (rank == rootrank) {
10359566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(sec0, &chart[0], &chart[1]));
1036c4762a1bSJed Brown   }
10379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(chart, 2, MPIU_INT, rootrank, comm));
10389566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(chart[1]-chart[0], &dofarr));
1039c4762a1bSJed Brown   if (rank == rootrank) {
1040c4762a1bSJed Brown     for (p = chart[0]; p < chart[1]; p++) {
10419566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sec0, p, &dofarr[p-chart[0]]));
1042c4762a1bSJed Brown     }
1043c4762a1bSJed Brown   }
10449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(dofarr, chart[1]-chart[0], MPIU_INT, rootrank, comm));
10459566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sec));
10469566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(sec, chart[0], chart[1]));
1047c4762a1bSJed Brown   for (p = chart[0]; p < chart[1]; p++) {
10489566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(sec, p, dofarr[p-chart[0]]));
1049c4762a1bSJed Brown   }
10509566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sec));
10519566063dSJacob Faibussowitsch   PetscCall(PetscFree(dofarr));
1052c4762a1bSJed Brown   *secout = sec;
1053c4762a1bSJed Brown   PetscFunctionReturn(0);
1054c4762a1bSJed Brown }
1055c4762a1bSJed Brown 
1056c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedVerticesCoordinatesToFaces_Private(DM ipdm, PortableBoundary bnd, IS *face_is)
1057c4762a1bSJed Brown {
1058c4762a1bSJed Brown   IS                  faces_expanded_is;
1059c4762a1bSJed Brown 
1060c4762a1bSJed Brown   PetscFunctionBegin;
10619566063dSJacob Faibussowitsch   PetscCall(DMPlexFindVertices(ipdm, bnd->coordinates, 0.0, &faces_expanded_is));
10629566063dSJacob Faibussowitsch   PetscCall(DMPlexExpandedVerticesToFaces_Private(ipdm, faces_expanded_is, bnd->depth, bnd->sections, face_is));
10639566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&faces_expanded_is));
1064c4762a1bSJed Brown   PetscFunctionReturn(0);
1065c4762a1bSJed Brown }
1066c4762a1bSJed Brown 
1067c4762a1bSJed Brown /* hack disabling DMPlexOrientInterface() call in DMPlexInterpolate() via -dm_plex_interpolate_orient_interfaces option */
1068c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM dm, PetscBool enable)
1069c4762a1bSJed Brown {
1070c4762a1bSJed Brown   PetscOptions      options = NULL;
1071c4762a1bSJed Brown   const char        *prefix = NULL;
1072c4762a1bSJed Brown   const char        opt[] = "-dm_plex_interpolate_orient_interfaces";
1073c4762a1bSJed Brown   char              prefix_opt[512];
1074c4762a1bSJed Brown   PetscBool         flg, set;
1075c4762a1bSJed Brown   static PetscBool  wasSetTrue = PETSC_FALSE;
1076c4762a1bSJed Brown 
1077c4762a1bSJed Brown   PetscFunctionBegin;
1078c4762a1bSJed Brown   if (dm) {
10799566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix));
1080c4762a1bSJed Brown     options = ((PetscObject)dm)->options;
1081c4762a1bSJed Brown   }
10829566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(prefix_opt, "-"));
10839566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(prefix_opt, prefix, sizeof(prefix_opt)));
10849566063dSJacob Faibussowitsch   PetscCall(PetscStrlcat(prefix_opt, &opt[1], sizeof(prefix_opt)));
10859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set));
1086c4762a1bSJed Brown   if (!enable) {
1087c4762a1bSJed Brown     if (set && flg) wasSetTrue = PETSC_TRUE;
10889566063dSJacob Faibussowitsch     PetscCall(PetscOptionsSetValue(options, prefix_opt, "0"));
1089c4762a1bSJed Brown   } else if (set && !flg) {
1090c4762a1bSJed Brown     if (wasSetTrue) {
10919566063dSJacob Faibussowitsch       PetscCall(PetscOptionsSetValue(options, prefix_opt, "1"));
1092c4762a1bSJed Brown     } else {
1093c4762a1bSJed Brown       /* default is PETSC_TRUE */
10949566063dSJacob Faibussowitsch       PetscCall(PetscOptionsClearValue(options, prefix_opt));
1095c4762a1bSJed Brown     }
1096c4762a1bSJed Brown     wasSetTrue = PETSC_FALSE;
1097c4762a1bSJed Brown   }
109876bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
10999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsGetBool(options, prefix, opt, &flg, &set));
1100*1dca8a05SBarry Smith     PetscCheck(!set || flg == enable,PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "PetscOptionsSetValue did not have the desired effect");
1101c4762a1bSJed Brown   }
1102c4762a1bSJed Brown   PetscFunctionReturn(0);
1103c4762a1bSJed Brown }
1104c4762a1bSJed Brown 
1105c4762a1bSJed Brown /* get coordinate description of the whole-domain boundary */
1106c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM dm, PortableBoundary *boundary)
1107c4762a1bSJed Brown {
1108c4762a1bSJed Brown   PortableBoundary       bnd0, bnd;
1109c4762a1bSJed Brown   MPI_Comm               comm;
1110c4762a1bSJed Brown   DM                     idm;
1111c4762a1bSJed Brown   DMLabel                label;
1112c4762a1bSJed Brown   PetscInt               d;
1113c4762a1bSJed Brown   const char             boundaryName[] = "DMPlexDistributeInterpolateMarkInterface_boundary";
1114c4762a1bSJed Brown   IS                     boundary_is;
1115c4762a1bSJed Brown   IS                     *boundary_expanded_iss;
1116c4762a1bSJed Brown   PetscMPIInt            rootrank = 0;
1117c4762a1bSJed Brown   PetscMPIInt            rank, size;
1118c4762a1bSJed Brown   PetscInt               value = 1;
1119c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1120c4762a1bSJed Brown   PetscBool              flg;
1121c4762a1bSJed Brown 
1122c4762a1bSJed Brown   PetscFunctionBegin;
11239566063dSJacob Faibussowitsch   PetscCall(PetscNew(&bnd));
11249566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
11259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
11269566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
11279566063dSJacob Faibussowitsch   PetscCall(DMPlexIsDistributed(dm, &flg));
112828b400f6SJacob Faibussowitsch   PetscCheck(!flg,comm, PETSC_ERR_ARG_WRONG, "serial DM (all points on one rank) needed");
1129c4762a1bSJed Brown 
1130c4762a1bSJed Brown   /* interpolate serial DM if not yet interpolated */
11319566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(dm, &intp));
1132c4762a1bSJed Brown   if (intp == DMPLEX_INTERPOLATED_FULL) {
1133c4762a1bSJed Brown     idm = dm;
11349566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)dm));
1135c4762a1bSJed Brown   } else {
11369566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(dm, &idm));
11379566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-idm_view"));
1138c4762a1bSJed Brown   }
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown   /* mark whole-domain boundary of the serial DM */
11419566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, boundaryName, &label));
11429566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(idm, label));
11439566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(idm, value, label));
11449566063dSJacob Faibussowitsch   PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-idm_boundary_view", comm));
11459566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumIS(label, value, &boundary_is));
1146c4762a1bSJed Brown 
1147c4762a1bSJed Brown   /* translate to coordinates */
11489566063dSJacob Faibussowitsch   PetscCall(PetscNew(&bnd0));
11499566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(idm));
1150c4762a1bSJed Brown   if (rank == rootrank) {
11519566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections));
11529566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocalTuple(dm, boundary_expanded_iss[0], NULL, &bnd0->coordinates));
1153c4762a1bSJed Brown     /* self-check */
1154c4762a1bSJed Brown     {
1155c4762a1bSJed Brown       IS is0;
11569566063dSJacob Faibussowitsch       PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(idm, bnd0, &is0));
11579566063dSJacob Faibussowitsch       PetscCall(ISEqual(is0, boundary_is, &flg));
115828b400f6SJacob Faibussowitsch       PetscCheck(flg,PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMPlexExpandedVerticesCoordinatesToFaces_Private produced a wrong IS");
11599566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is0));
1160c4762a1bSJed Brown     }
1161c4762a1bSJed Brown   } else {
11629566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &bnd0->coordinates));
1163c4762a1bSJed Brown   }
1164c4762a1bSJed Brown 
1165c4762a1bSJed Brown   {
1166c4762a1bSJed Brown     Vec         tmp;
1167c4762a1bSJed Brown     VecScatter  sc;
1168c4762a1bSJed Brown     IS          xis;
1169c4762a1bSJed Brown     PetscInt    n;
1170c4762a1bSJed Brown 
1171c4762a1bSJed Brown     /* just convert seq vectors to mpi vector */
11729566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(bnd0->coordinates, &n));
11739566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&n, 1, MPIU_INT, rootrank, comm));
1174c4762a1bSJed Brown     if (rank == rootrank) {
11759566063dSJacob Faibussowitsch       PetscCall(VecCreateMPI(comm, n, n, &tmp));
1176c4762a1bSJed Brown     } else {
11779566063dSJacob Faibussowitsch       PetscCall(VecCreateMPI(comm, 0, n, &tmp));
1178c4762a1bSJed Brown     }
11799566063dSJacob Faibussowitsch     PetscCall(VecCopy(bnd0->coordinates, tmp));
11809566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&bnd0->coordinates));
1181c4762a1bSJed Brown     bnd0->coordinates = tmp;
1182c4762a1bSJed Brown 
1183c4762a1bSJed Brown     /* replicate coordinates from root rank to all ranks */
11849566063dSJacob Faibussowitsch     PetscCall(VecCreateMPI(comm, n, n*size, &bnd->coordinates));
11859566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, n, 0, 1, &xis));
11869566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(bnd0->coordinates, xis, bnd->coordinates, NULL, &sc));
11879566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD));
11889566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(  sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD));
11899566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&sc));
11909566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&xis));
1191c4762a1bSJed Brown   }
1192c4762a1bSJed Brown   bnd->depth = bnd0->depth;
11939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&bnd->depth, 1, MPIU_INT, rootrank, comm));
11949566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(bnd->depth, &bnd->sections));
1195c4762a1bSJed Brown   for (d=0; d<bnd->depth; d++) {
11969566063dSJacob Faibussowitsch     PetscCall(PetscSectionReplicate_Private(comm, rootrank, (rank == rootrank) ? bnd0->sections[d] : NULL, &bnd->sections[d]));
1197c4762a1bSJed Brown   }
1198c4762a1bSJed Brown 
1199c4762a1bSJed Brown   if (rank == rootrank) {
12009566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections));
1201c4762a1bSJed Brown   }
12029566063dSJacob Faibussowitsch   PetscCall(PortableBoundaryDestroy(&bnd0));
12039566063dSJacob Faibussowitsch   PetscCall(DMRemoveLabelBySelf(idm, &label, PETSC_TRUE));
12049566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
12059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&boundary_is));
12069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&idm));
1207c4762a1bSJed Brown   *boundary = bnd;
1208c4762a1bSJed Brown   PetscFunctionReturn(0);
1209c4762a1bSJed Brown }
1210c4762a1bSJed Brown 
1211c4762a1bSJed Brown /* get faces of inter-partition interface */
1212c4762a1bSJed Brown static PetscErrorCode DMPlexGetInterfaceFaces_Private(DM ipdm, IS boundary_faces_is, IS *interface_faces_is)
1213c4762a1bSJed Brown {
1214c4762a1bSJed Brown   MPI_Comm               comm;
1215c4762a1bSJed Brown   DMLabel                label;
1216c4762a1bSJed Brown   IS                     part_boundary_faces_is;
1217c4762a1bSJed Brown   const char             partBoundaryName[] = "DMPlexDistributeInterpolateMarkInterface_partBoundary";
1218c4762a1bSJed Brown   PetscInt               value = 1;
1219c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1220c4762a1bSJed Brown 
1221c4762a1bSJed Brown   PetscFunctionBegin;
12229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm));
12239566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp));
122408401ef6SPierre Jolivet   PetscCheck(intp == DMPLEX_INTERPOLATED_FULL,comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex");
1225c4762a1bSJed Brown 
1226c4762a1bSJed Brown   /* get ipdm partition boundary (partBoundary) */
12279566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, partBoundaryName, &label));
12289566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(ipdm, label));
12299566063dSJacob Faibussowitsch   PetscCall(DMPlexMarkBoundaryFaces(ipdm, value, label));
12309566063dSJacob Faibussowitsch   PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-ipdm_part_boundary_view", comm));
12319566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, &part_boundary_faces_is));
12329566063dSJacob Faibussowitsch   PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE));
12339566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
1234c4762a1bSJed Brown 
1235c4762a1bSJed Brown   /* remove ipdm whole-domain boundary (boundary_faces_is) from ipdm partition boundary (part_boundary_faces_is), resulting just in inter-partition interface */
12369566063dSJacob Faibussowitsch   PetscCall(ISDifference(part_boundary_faces_is,boundary_faces_is,interface_faces_is));
12379566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&part_boundary_faces_is));
1238c4762a1bSJed Brown   PetscFunctionReturn(0);
1239c4762a1bSJed Brown }
1240c4762a1bSJed Brown 
1241c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices */
1242c4762a1bSJed Brown static PetscErrorCode DMPlexComputeCompleteInterface_Private(DM ipdm, IS interface_faces_is, IS *interface_is)
1243c4762a1bSJed Brown {
1244c4762a1bSJed Brown   DMLabel                label;
1245c4762a1bSJed Brown   PetscInt               value = 1;
1246c4762a1bSJed Brown   const char             interfaceName[] = "DMPlexDistributeInterpolateMarkInterface_interface";
1247c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1248c4762a1bSJed Brown   MPI_Comm               comm;
1249c4762a1bSJed Brown 
1250c4762a1bSJed Brown   PetscFunctionBegin;
12519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm));
12529566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(ipdm, &intp));
125308401ef6SPierre Jolivet   PetscCheck(intp == DMPLEX_INTERPOLATED_FULL,comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex");
1254c4762a1bSJed Brown 
12559566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, interfaceName, &label));
12569566063dSJacob Faibussowitsch   PetscCall(DMAddLabel(ipdm, label));
12579566063dSJacob Faibussowitsch   PetscCall(DMLabelSetStratumIS(label, value, interface_faces_is));
12589566063dSJacob Faibussowitsch   PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_faces_view", comm));
12599566063dSJacob Faibussowitsch   PetscCall(DMPlexLabelComplete(ipdm, label));
12609566063dSJacob Faibussowitsch   PetscCall(DMLabelViewFromOptionsOnComm_Private(label, "-interface_view", comm));
12619566063dSJacob Faibussowitsch   PetscCall(DMLabelGetStratumISOnComm_Private(label, value, comm, interface_is));
12629566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*interface_is, "interface_is"));
12639566063dSJacob Faibussowitsch   PetscCall(ISViewFromOptions(*interface_is, NULL, "-interface_is_view"));
12649566063dSJacob Faibussowitsch   PetscCall(DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE));
12659566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&label));
1266c4762a1bSJed Brown   PetscFunctionReturn(0);
1267c4762a1bSJed Brown }
1268c4762a1bSJed Brown 
1269c4762a1bSJed Brown static PetscErrorCode PointSFGetOutwardInterfacePoints(PetscSF sf, IS *is)
1270c4762a1bSJed Brown {
1271c4762a1bSJed Brown   PetscInt        n;
1272c4762a1bSJed Brown   const PetscInt  *arr;
1273c4762a1bSJed Brown 
1274c4762a1bSJed Brown   PetscFunctionBegin;
12759566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, NULL, &n, &arr, NULL));
12769566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_USE_POINTER, is));
1277c4762a1bSJed Brown   PetscFunctionReturn(0);
1278c4762a1bSJed Brown }
1279c4762a1bSJed Brown 
1280c4762a1bSJed Brown static PetscErrorCode PointSFGetInwardInterfacePoints(PetscSF sf, IS *is)
1281c4762a1bSJed Brown {
1282c4762a1bSJed Brown   PetscInt        n;
1283c4762a1bSJed Brown   const PetscInt  *rootdegree;
1284c4762a1bSJed Brown   PetscInt        *arr;
1285c4762a1bSJed Brown 
1286c4762a1bSJed Brown   PetscFunctionBegin;
12879566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sf));
12889566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(sf, &rootdegree));
12899566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(sf, &rootdegree));
12909566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeMultiRootOriginalNumbering(sf, rootdegree, &n, &arr));
12919566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_OWN_POINTER, is));
1292c4762a1bSJed Brown   PetscFunctionReturn(0);
1293c4762a1bSJed Brown }
1294c4762a1bSJed Brown 
1295c4762a1bSJed Brown static PetscErrorCode PointSFGetInterfacePoints_Private(PetscSF pointSF, IS *is)
1296c4762a1bSJed Brown {
1297c4762a1bSJed Brown   IS pointSF_out_is, pointSF_in_is;
1298c4762a1bSJed Brown 
1299c4762a1bSJed Brown   PetscFunctionBegin;
13009566063dSJacob Faibussowitsch   PetscCall(PointSFGetOutwardInterfacePoints(pointSF, &pointSF_out_is));
13019566063dSJacob Faibussowitsch   PetscCall(PointSFGetInwardInterfacePoints(pointSF, &pointSF_in_is));
13029566063dSJacob Faibussowitsch   PetscCall(ISExpand(pointSF_out_is, pointSF_in_is, is));
13039566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointSF_out_is));
13049566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointSF_in_is));
1305c4762a1bSJed Brown   PetscFunctionReturn(0);
1306c4762a1bSJed Brown }
1307c4762a1bSJed Brown 
13085f80ce2aSJacob Faibussowitsch #define CHKERRMY(ierr) PetscCheck(!ierr,PETSC_COMM_SELF, PETSC_ERR_PLIB, "PointSF is wrong. Unable to show details!")
1309c4762a1bSJed Brown 
1310c4762a1bSJed Brown static PetscErrorCode ViewPointsWithType_Internal(DM dm, IS pointsIS, PetscViewer v)
1311c4762a1bSJed Brown {
1312c4762a1bSJed Brown   DMLabel         label;
1313c4762a1bSJed Brown   PetscSection    coordsSection;
1314c4762a1bSJed Brown   Vec             coordsVec;
1315c4762a1bSJed Brown   PetscScalar     *coordsScalar;
1316c4762a1bSJed Brown   PetscInt        coneSize, depth, dim, i, p, npoints;
1317c4762a1bSJed Brown   const PetscInt  *points;
1318c4762a1bSJed Brown 
1319c4762a1bSJed Brown   PetscFunctionBegin;
13209566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
13219566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &coordsSection));
13229566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &coordsVec));
13239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordsVec, &coordsScalar));
13249566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(pointsIS, &npoints));
13259566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(pointsIS, &points));
13269566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &label));
13279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(v));
1328c4762a1bSJed Brown   for (i=0; i<npoints; i++) {
1329c4762a1bSJed Brown     p = points[i];
13309566063dSJacob Faibussowitsch     PetscCall(DMLabelGetValue(label, p, &depth));
1331c4762a1bSJed Brown     if (!depth) {
1332d3e1f4ccSVaclav Hapla       PetscInt        n, o;
1333c4762a1bSJed Brown       char            coordstr[128];
1334c4762a1bSJed Brown 
13359566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(coordsSection, p, &n));
13369566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordsSection, p, &o));
13379566063dSJacob Faibussowitsch       PetscCall(coord2str(coordstr, sizeof(coordstr), n, &coordsScalar[o], 1.0));
133863a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(v, "vertex %" PetscInt_FMT " w/ coordinates %s\n", p, coordstr));
1339c4762a1bSJed Brown     } else {
1340c4762a1bSJed Brown       char            entityType[16];
1341c4762a1bSJed Brown 
1342c4762a1bSJed Brown       switch (depth) {
13439566063dSJacob Faibussowitsch         case 1: PetscCall(PetscStrcpy(entityType, "edge")); break;
13449566063dSJacob Faibussowitsch         case 2: PetscCall(PetscStrcpy(entityType, "face")); break;
13459566063dSJacob Faibussowitsch         case 3: PetscCall(PetscStrcpy(entityType, "cell")); break;
13462a27bf02SStefano Zampini         default: SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for depth <= 3");
1347c4762a1bSJed Brown       }
1348c4762a1bSJed Brown       if (depth == dim && dim < 3) {
13499566063dSJacob Faibussowitsch         PetscCall(PetscStrlcat(entityType, " (cell)", sizeof(entityType)));
1350c4762a1bSJed Brown       }
135163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(v, "%s %" PetscInt_FMT "\n", entityType, p));
1352c4762a1bSJed Brown     }
13539566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
1354c4762a1bSJed Brown     if (coneSize) {
1355c4762a1bSJed Brown       const PetscInt *cone;
1356c4762a1bSJed Brown       IS             coneIS;
1357c4762a1bSJed Brown 
13589566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCone(dm, p, &cone));
13599566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, coneSize, cone, PETSC_USE_POINTER, &coneIS));
13609566063dSJacob Faibussowitsch       PetscCall(ViewPointsWithType_Internal(dm, coneIS, v));
13619566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&coneIS));
1362c4762a1bSJed Brown     }
1363c4762a1bSJed Brown   }
13649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(v));
13659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordsVec, &coordsScalar));
13669566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(pointsIS, &points));
1367c4762a1bSJed Brown   PetscFunctionReturn(0);
1368c4762a1bSJed Brown }
1369c4762a1bSJed Brown 
1370c4762a1bSJed Brown static PetscErrorCode ViewPointsWithType(DM dm, IS points, PetscViewer v)
1371c4762a1bSJed Brown {
1372c4762a1bSJed Brown   PetscBool       flg;
1373c4762a1bSJed Brown   PetscInt        npoints;
1374c4762a1bSJed Brown   PetscMPIInt     rank;
1375c4762a1bSJed Brown 
1376c4762a1bSJed Brown   PetscFunctionBegin;
13779566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &flg));
137828b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for ASCII viewer");
13799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)v), &rank));
13809566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushSynchronized(v));
13819566063dSJacob Faibussowitsch   PetscCall(ISGetLocalSize(points, &npoints));
1382c4762a1bSJed Brown   if (npoints) {
13839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", rank));
13849566063dSJacob Faibussowitsch     PetscCall(ViewPointsWithType_Internal(dm, points, v));
1385c4762a1bSJed Brown   }
13869566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(v));
13879566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopSynchronized(v));
1388c4762a1bSJed Brown   PetscFunctionReturn(0);
1389c4762a1bSJed Brown }
1390c4762a1bSJed Brown 
1391c4762a1bSJed Brown static PetscErrorCode DMPlexComparePointSFWithInterface_Private(DM ipdm, IS interface_is)
1392c4762a1bSJed Brown {
1393c4762a1bSJed Brown   PetscSF         pointsf;
1394c4762a1bSJed Brown   IS              pointsf_is;
1395c4762a1bSJed Brown   PetscBool       flg;
1396c4762a1bSJed Brown   MPI_Comm        comm;
1397c4762a1bSJed Brown   PetscMPIInt     size;
1398c4762a1bSJed Brown 
1399c4762a1bSJed Brown   PetscFunctionBegin;
14009566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ipdm, &comm));
14019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
14029566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(ipdm, &pointsf));
1403c4762a1bSJed Brown   if (pointsf) {
1404c4762a1bSJed Brown     PetscInt nroots;
14059566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(pointsf, &nroots, NULL, NULL, NULL));
1406c4762a1bSJed Brown     if (nroots < 0) pointsf = NULL; /* uninitialized SF */
1407c4762a1bSJed Brown   }
1408c4762a1bSJed Brown   if (!pointsf) {
1409c4762a1bSJed Brown     PetscInt N=0;
14109566063dSJacob Faibussowitsch     if (interface_is) PetscCall(ISGetSize(interface_is, &N));
141128b400f6SJacob Faibussowitsch     PetscCheck(!N,comm, PETSC_ERR_PLIB, "interface_is should be NULL or empty for PointSF being NULL");
1412c4762a1bSJed Brown     PetscFunctionReturn(0);
1413c4762a1bSJed Brown   }
1414c4762a1bSJed Brown 
1415c4762a1bSJed Brown   /* get PointSF points as IS pointsf_is */
14169566063dSJacob Faibussowitsch   PetscCall(PointSFGetInterfacePoints_Private(pointsf, &pointsf_is));
1417c4762a1bSJed Brown 
1418c4762a1bSJed Brown   /* compare pointsf_is with interface_is */
14199566063dSJacob Faibussowitsch   PetscCall(ISEqual(interface_is, pointsf_is, &flg));
14209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE,&flg,1,MPIU_BOOL,MPI_LAND,comm));
1421c4762a1bSJed Brown   if (!flg) {
1422c4762a1bSJed Brown     IS pointsf_extra_is, pointsf_missing_is;
1423c4762a1bSJed Brown     PetscViewer errv = PETSC_VIEWER_STDERR_(comm);
14245f80ce2aSJacob Faibussowitsch     CHKERRMY(ISDifference(interface_is, pointsf_is, &pointsf_missing_is));
14255f80ce2aSJacob Faibussowitsch     CHKERRMY(ISDifference(pointsf_is, interface_is, &pointsf_extra_is));
14265f80ce2aSJacob Faibussowitsch     CHKERRMY(PetscViewerASCIIPrintf(errv, "Points missing in PointSF:\n"));
14275f80ce2aSJacob Faibussowitsch     CHKERRMY(ViewPointsWithType(ipdm, pointsf_missing_is, errv));
14285f80ce2aSJacob Faibussowitsch     CHKERRMY(PetscViewerASCIIPrintf(errv, "Extra points in PointSF:\n"));
14295f80ce2aSJacob Faibussowitsch     CHKERRMY(ViewPointsWithType(ipdm, pointsf_extra_is, errv));
14305f80ce2aSJacob Faibussowitsch     CHKERRMY(ISDestroy(&pointsf_extra_is));
14315f80ce2aSJacob Faibussowitsch     CHKERRMY(ISDestroy(&pointsf_missing_is));
1432c4762a1bSJed Brown     SETERRQ(comm, PETSC_ERR_PLIB, "PointSF is wrong! See details above.");
1433c4762a1bSJed Brown   }
14349566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&pointsf_is));
1435c4762a1bSJed Brown   PetscFunctionReturn(0);
1436c4762a1bSJed Brown }
1437c4762a1bSJed Brown 
1438c4762a1bSJed Brown /* remove faces & edges from label, leave just vertices */
1439c4762a1bSJed Brown static PetscErrorCode DMPlexISFilterVertices_Private(DM dm, IS points)
1440c4762a1bSJed Brown {
1441c4762a1bSJed Brown   PetscInt        vStart, vEnd;
1442c4762a1bSJed Brown   MPI_Comm        comm;
1443c4762a1bSJed Brown 
1444c4762a1bSJed Brown   PetscFunctionBegin;
14459566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
14469566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
14479566063dSJacob Faibussowitsch   PetscCall(ISGeneralFilter(points, vStart, vEnd));
1448c4762a1bSJed Brown   PetscFunctionReturn(0);
1449c4762a1bSJed Brown }
1450c4762a1bSJed Brown 
1451c4762a1bSJed Brown /*
14522e58f0efSBarry Smith   DMPlexCheckPointSFHeavy - Thoroughly test that the PointSF after parallel DMPlexInterpolate() includes exactly all interface points.
1453c4762a1bSJed Brown 
1454c4762a1bSJed Brown   Collective
1455c4762a1bSJed Brown 
1456c4762a1bSJed Brown   Input Parameters:
1457c4762a1bSJed Brown . dm - The DMPlex object
1458c4762a1bSJed Brown 
1459c4762a1bSJed Brown   Notes:
1460c4762a1bSJed Brown   The input DMPlex must be serial (one partition has all points, the other partitions have no points).
1461c4762a1bSJed Brown   This is a heavy test which involves DMPlexInterpolate() if the input DM is not interpolated yet, and depends on having a representation of the whole-domain boundary (PortableBoundary), which can be obtained only by DMPlexGetExpandedBoundary_Private() (which involves DMPlexInterpolate() of a sequential DM).
1462c4762a1bSJed Brown   This is mainly intended for debugging/testing purposes.
1463c4762a1bSJed Brown 
1464c4762a1bSJed Brown   Algorithm:
1465c4762a1bSJed Brown   1. boundary faces of the serial version of the whole mesh are found using DMPlexMarkBoundaryFaces()
1466c4762a1bSJed Brown   2. boundary faces are translated into vertices using DMPlexGetConeRecursive() and these are translated into coordinates - this description (aka PortableBoundary) is completely independent of partitioning and point numbering
1467c4762a1bSJed Brown   3. the mesh is distributed or loaded in parallel
1468c4762a1bSJed Brown   4. boundary faces of the distributed mesh are reconstructed from PortableBoundary using DMPlexFindVertices()
1469c4762a1bSJed Brown   5. partition boundary faces of the parallel mesh are found using DMPlexMarkBoundaryFaces()
1470c4762a1bSJed Brown   6. partition interfaces are computed as set difference of partition boundary faces minus the reconstructed boundary
1471c4762a1bSJed Brown   7. check that interface covered by PointSF (union of inward and outward points) is equal to the partition interface for each rank, otherwise print the difference and throw an error
1472c4762a1bSJed Brown 
1473c4762a1bSJed Brown   Level: developer
1474c4762a1bSJed Brown 
1475c4762a1bSJed Brown .seealso: DMGetPointSF(), DMPlexCheckSymmetry(), DMPlexCheckSkeleton(), DMPlexCheckFaces()
1476c4762a1bSJed Brown */
1477c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM dm, PortableBoundary bnd)
1478c4762a1bSJed Brown {
1479c4762a1bSJed Brown   DM                     ipdm=NULL;
1480c4762a1bSJed Brown   IS                     boundary_faces_is, interface_faces_is, interface_is;
1481c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1482c4762a1bSJed Brown   MPI_Comm               comm;
1483c4762a1bSJed Brown 
1484c4762a1bSJed Brown   PetscFunctionBegin;
14859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1486c4762a1bSJed Brown 
14879566063dSJacob Faibussowitsch   PetscCall(DMPlexIsInterpolatedCollective(dm, &intp));
1488c4762a1bSJed Brown   if (intp == DMPLEX_INTERPOLATED_FULL) {
1489c4762a1bSJed Brown     ipdm = dm;
1490c4762a1bSJed Brown   } else {
1491c4762a1bSJed Brown     /* create temporary interpolated DM if input DM is not interpolated */
14929566063dSJacob Faibussowitsch     PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_FALSE));
14939566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(dm, &ipdm)); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexComparePointSFWithInterface_Private() below */
14949566063dSJacob Faibussowitsch     PetscCall(DMPlexSetOrientInterface_Private(dm, PETSC_TRUE));
1495c4762a1bSJed Brown   }
14969566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(ipdm, NULL, "-ipdm_view"));
1497c4762a1bSJed Brown 
1498c4762a1bSJed Brown   /* recover ipdm whole-domain boundary faces from the expanded vertices coordinates */
14999566063dSJacob Faibussowitsch   PetscCall(DMPlexExpandedVerticesCoordinatesToFaces_Private(ipdm, bnd, &boundary_faces_is));
1500c4762a1bSJed Brown   /* get inter-partition interface faces (interface_faces_is)*/
15019566063dSJacob Faibussowitsch   PetscCall(DMPlexGetInterfaceFaces_Private(ipdm, boundary_faces_is, &interface_faces_is));
1502c4762a1bSJed Brown   /* compute inter-partition interface including edges and vertices (interface_is) */
15039566063dSJacob Faibussowitsch   PetscCall(DMPlexComputeCompleteInterface_Private(ipdm, interface_faces_is, &interface_is));
1504c4762a1bSJed Brown   /* destroy immediate ISs */
15059566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&boundary_faces_is));
15069566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&interface_faces_is));
1507c4762a1bSJed Brown 
1508c4762a1bSJed Brown   /* for uninterpolated case, keep just vertices in interface */
1509c4762a1bSJed Brown   if (!intp) {
15109566063dSJacob Faibussowitsch     PetscCall(DMPlexISFilterVertices_Private(ipdm, interface_is));
15119566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&ipdm));
1512c4762a1bSJed Brown   }
1513c4762a1bSJed Brown 
1514c4762a1bSJed Brown   /* compare PointSF with the boundary reconstructed from coordinates */
15159566063dSJacob Faibussowitsch   PetscCall(DMPlexComparePointSFWithInterface_Private(dm, interface_is));
15169566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "DMPlexCheckPointSFHeavy PASSED\n"));
15179566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&interface_is));
1518c4762a1bSJed Brown   PetscFunctionReturn(0);
1519c4762a1bSJed Brown }
1520c4762a1bSJed Brown 
1521c4762a1bSJed Brown int main(int argc, char **argv)
1522c4762a1bSJed Brown {
1523c4762a1bSJed Brown   DM             dm;
1524c4762a1bSJed Brown   AppCtx         user;
1525c4762a1bSJed Brown 
15269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15279566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("create",&stage[0]));
15289566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("distribute",&stage[1]));
15299566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("interpolate",&stage[2]));
15309566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
15319566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
1532c4762a1bSJed Brown   if (user.nPointsToExpand) {
15339566063dSJacob Faibussowitsch     PetscCall(TestExpandPoints(dm, &user));
1534c4762a1bSJed Brown   }
1535c4762a1bSJed Brown   if (user.ncoords) {
1536d3e1f4ccSVaclav Hapla     Vec coords;
1537d3e1f4ccSVaclav Hapla 
15389566063dSJacob Faibussowitsch     PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, user.ncoords, user.ncoords, user.coords, &coords));
15399566063dSJacob Faibussowitsch     PetscCall(ViewVerticesFromCoords(dm, coords, user.coordsTol, PETSC_VIEWER_STDOUT_WORLD));
15409566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&coords));
1541c4762a1bSJed Brown   }
15429566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
15439566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1544b122ec5aSJacob Faibussowitsch   return 0;
1545c4762a1bSJed Brown }
1546c4762a1bSJed Brown 
1547c4762a1bSJed Brown /*TEST
1548c4762a1bSJed Brown 
1549c4762a1bSJed Brown   testset:
1550c4762a1bSJed Brown     nsize: 2
1551c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail
1552c4762a1bSJed Brown     args: -dm_plex_check_all
1553c4762a1bSJed Brown     test:
1554c4762a1bSJed Brown       suffix: 1_tri_dist0
1555c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1556c4762a1bSJed Brown     test:
1557c4762a1bSJed Brown       suffix: 1_tri_dist1
1558c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1559c4762a1bSJed Brown     test:
1560c4762a1bSJed Brown       suffix: 1_quad_dist0
1561c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output}
1562c4762a1bSJed Brown     test:
1563c4762a1bSJed Brown       suffix: 1_quad_dist1
1564c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output}
1565c4762a1bSJed Brown     test:
1566c4762a1bSJed Brown       suffix: 1_1d_dist1
1567c4762a1bSJed Brown       args: -dim 1 -distribute 1
1568c4762a1bSJed Brown 
1569c4762a1bSJed Brown   testset:
1570c4762a1bSJed Brown     nsize: 3
1571c4762a1bSJed Brown     args: -testnum 1 -interpolate create
1572c4762a1bSJed Brown     args: -dm_plex_check_all
1573c4762a1bSJed Brown     test:
1574c4762a1bSJed Brown       suffix: 2
1575c4762a1bSJed Brown       args: -dm_view ascii::ascii_info_detail
1576c4762a1bSJed Brown     test:
1577c4762a1bSJed Brown       suffix: 2a
1578c4762a1bSJed Brown       args: -dm_plex_check_cones_conform_on_interfaces_verbose
1579c4762a1bSJed Brown     test:
1580c4762a1bSJed Brown       suffix: 2b
1581c4762a1bSJed Brown       args: -test_expand_points 0,1,2,5,6
1582c4762a1bSJed Brown     test:
1583c4762a1bSJed Brown       suffix: 2c
1584c4762a1bSJed Brown       args: -test_expand_points 0,1,2,5,6 -test_expand_points_empty
1585c4762a1bSJed Brown 
1586c4762a1bSJed Brown   testset:
1587c4762a1bSJed Brown     # the same as 1% for 3D
1588c4762a1bSJed Brown     nsize: 2
1589c4762a1bSJed Brown     args: -dim 3 -dm_view ascii::ascii_info_detail
1590c4762a1bSJed Brown     args: -dm_plex_check_all
1591c4762a1bSJed Brown     test:
1592c4762a1bSJed Brown       suffix: 4_tet_dist0
1593c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1594c4762a1bSJed Brown     test:
1595c4762a1bSJed Brown       suffix: 4_tet_dist1
1596c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1597c4762a1bSJed Brown     test:
1598c4762a1bSJed Brown       suffix: 4_hex_dist0
1599c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output}
1600c4762a1bSJed Brown     test:
1601c4762a1bSJed Brown       suffix: 4_hex_dist1
1602c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output}
1603c4762a1bSJed Brown 
1604c4762a1bSJed Brown   test:
1605c4762a1bSJed Brown     # the same as 4_tet_dist0 but test different initial orientations
1606c4762a1bSJed Brown     suffix: 4_tet_test_orient
1607c4762a1bSJed Brown     nsize: 2
1608c4762a1bSJed Brown     args: -dim 3 -distribute 0
1609c4762a1bSJed Brown     args: -dm_plex_check_all
1610c4762a1bSJed Brown     args: -rotate_interface_0 {{0 1 2 11 12 13}}
1611c4762a1bSJed Brown     args: -rotate_interface_1 {{0 1 2 11 12 13}}
1612c4762a1bSJed Brown 
1613c4762a1bSJed Brown   testset:
1614c4762a1bSJed Brown     requires: exodusii
1615c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo
1616c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail
1617c4762a1bSJed Brown     args: -dm_plex_check_all
1618c4762a1bSJed Brown     args: -custom_view
1619c4762a1bSJed Brown     test:
1620c4762a1bSJed Brown       suffix: 5_seq
1621c4762a1bSJed Brown       nsize: 1
1622c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1623c4762a1bSJed Brown     test:
162430602db0SMatthew G. Knepley       # Detail viewing in a non-distributed mesh is broken because the DMLabelView() is collective, but the label is not shared
1625c4762a1bSJed Brown       suffix: 5_dist0
1626c4762a1bSJed Brown       nsize: 2
162730602db0SMatthew G. Knepley       args: -distribute 0 -interpolate {{none create}separate output} -dm_view
1628c4762a1bSJed Brown     test:
1629c4762a1bSJed Brown       suffix: 5_dist1
1630c4762a1bSJed Brown       nsize: 2
1631c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1632c4762a1bSJed Brown 
1633c4762a1bSJed Brown   testset:
1634c4762a1bSJed Brown     nsize: {{1 2 4}}
1635c4762a1bSJed Brown     args: -use_generator
1636c4762a1bSJed Brown     args: -dm_plex_check_all
1637c4762a1bSJed Brown     args: -distribute -interpolate none
1638c4762a1bSJed Brown     test:
1639c4762a1bSJed Brown       suffix: 6_tri
1640c4762a1bSJed Brown       requires: triangle
1641c0517cd5SMatthew G. Knepley       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_generator triangle
1642c4762a1bSJed Brown     test:
1643c4762a1bSJed Brown       suffix: 6_quad
1644c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1645c4762a1bSJed Brown     test:
1646c4762a1bSJed Brown       suffix: 6_tet
1647c4762a1bSJed Brown       requires: ctetgen
1648c0517cd5SMatthew G. Knepley       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_generator ctetgen
1649c4762a1bSJed Brown     test:
1650c4762a1bSJed Brown       suffix: 6_hex
1651c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1652c4762a1bSJed Brown   testset:
1653c4762a1bSJed Brown     nsize: {{1 2 4}}
1654c4762a1bSJed Brown     args: -use_generator
1655c4762a1bSJed Brown     args: -dm_plex_check_all
1656c4762a1bSJed Brown     args: -distribute -interpolate create
1657c4762a1bSJed Brown     test:
1658c4762a1bSJed Brown       suffix: 6_int_tri
1659c4762a1bSJed Brown       requires: triangle
1660c0517cd5SMatthew G. Knepley       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_generator triangle
1661c4762a1bSJed Brown     test:
1662c4762a1bSJed Brown       suffix: 6_int_quad
1663c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1664c4762a1bSJed Brown     test:
1665c4762a1bSJed Brown       suffix: 6_int_tet
1666c4762a1bSJed Brown       requires: ctetgen
1667c0517cd5SMatthew G. Knepley       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_generator ctetgen
1668c4762a1bSJed Brown     test:
1669c4762a1bSJed Brown       suffix: 6_int_hex
1670c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1671c4762a1bSJed Brown   testset:
1672c4762a1bSJed Brown     nsize: {{2 4}}
1673c4762a1bSJed Brown     args: -use_generator
1674c4762a1bSJed Brown     args: -dm_plex_check_all
1675c4762a1bSJed Brown     args: -distribute -interpolate after_distribute
1676c4762a1bSJed Brown     test:
1677c4762a1bSJed Brown       suffix: 6_parint_tri
1678c4762a1bSJed Brown       requires: triangle
1679c0517cd5SMatthew G. Knepley       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_generator triangle
1680c4762a1bSJed Brown     test:
1681c4762a1bSJed Brown       suffix: 6_parint_quad
1682c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1683c4762a1bSJed Brown     test:
1684c4762a1bSJed Brown       suffix: 6_parint_tet
1685c4762a1bSJed Brown       requires: ctetgen
1686c0517cd5SMatthew G. Knepley       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_generator ctetgen
1687c4762a1bSJed Brown     test:
1688c4762a1bSJed Brown       suffix: 6_parint_hex
1689c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1690c4762a1bSJed Brown 
1691c4762a1bSJed Brown   testset: # 7 EXODUS
1692c4762a1bSJed Brown     requires: exodusii
1693c4762a1bSJed Brown     args: -dm_plex_check_all
1694c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
1695c4762a1bSJed Brown     args: -distribute
1696c4762a1bSJed Brown     test: # seq load, simple partitioner
1697c4762a1bSJed Brown       suffix: 7_exo
1698c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1699c4762a1bSJed Brown       args: -interpolate none
1700c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1701c4762a1bSJed Brown       suffix: 7_exo_int_simple
1702c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1703c4762a1bSJed Brown       args: -interpolate create
1704c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1705c4762a1bSJed Brown       suffix: 7_exo_int_metis
1706c4762a1bSJed Brown       requires: parmetis
1707c4762a1bSJed Brown       nsize: {{2 4 5}}
1708c4762a1bSJed Brown       args: -interpolate create
1709c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1710c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1711c4762a1bSJed Brown       suffix: 7_exo_simple_int
1712c4762a1bSJed Brown       nsize: {{2 4 5}}
1713c4762a1bSJed Brown       args: -interpolate after_distribute
1714c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1715c4762a1bSJed Brown       suffix: 7_exo_metis_int
1716c4762a1bSJed Brown       requires: parmetis
1717c4762a1bSJed Brown       nsize: {{2 4 5}}
1718c4762a1bSJed Brown       args: -interpolate after_distribute
1719c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1720c4762a1bSJed Brown 
1721c4762a1bSJed Brown   testset: # 7 HDF5 SEQUANTIAL LOAD
1722c4762a1bSJed Brown     requires: hdf5 !complex
1723c4762a1bSJed Brown     args: -dm_plex_check_all
1724c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1725c4762a1bSJed Brown     args: -dm_plex_hdf5_force_sequential
1726c4762a1bSJed Brown     args: -distribute
1727c4762a1bSJed Brown     test: # seq load, simple partitioner
1728c4762a1bSJed Brown       suffix: 7_seq_hdf5_simple
1729c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1730c4762a1bSJed Brown       args: -interpolate none
1731c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1732c4762a1bSJed Brown       suffix: 7_seq_hdf5_int_simple
1733c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1734c4762a1bSJed Brown       args: -interpolate after_create
1735c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1736c4762a1bSJed Brown       nsize: {{2 4 5}}
1737c4762a1bSJed Brown       suffix: 7_seq_hdf5_int_metis
1738c4762a1bSJed Brown       requires: parmetis
1739c4762a1bSJed Brown       args: -interpolate after_create
1740c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1741c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1742c4762a1bSJed Brown       suffix: 7_seq_hdf5_simple_int
1743c4762a1bSJed Brown       nsize: {{2 4 5}}
1744c4762a1bSJed Brown       args: -interpolate after_distribute
1745c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1746c4762a1bSJed Brown       nsize: {{2 4 5}}
1747c4762a1bSJed Brown       suffix: 7_seq_hdf5_metis_int
1748c4762a1bSJed Brown       requires: parmetis
1749c4762a1bSJed Brown       args: -interpolate after_distribute
1750c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1751c4762a1bSJed Brown 
1752c4762a1bSJed Brown   testset: # 7 HDF5 PARALLEL LOAD
1753c4762a1bSJed Brown     requires: hdf5 !complex
1754c4762a1bSJed Brown     nsize: {{2 4 5}}
1755c4762a1bSJed Brown     args: -dm_plex_check_all
1756c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1757c4762a1bSJed Brown     test: # par load
1758c4762a1bSJed Brown       suffix: 7_par_hdf5
1759c4762a1bSJed Brown       args: -interpolate none
1760c4762a1bSJed Brown     test: # par load, par interpolation
1761c4762a1bSJed Brown       suffix: 7_par_hdf5_int
1762c4762a1bSJed Brown       args: -interpolate after_create
1763c4762a1bSJed Brown     test: # par load, parmetis repartitioner
1764c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1765c4762a1bSJed Brown       suffix: 7_par_hdf5_parmetis
1766c4762a1bSJed Brown       requires: parmetis
1767c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1768c4762a1bSJed Brown       args: -interpolate none
1769c4762a1bSJed Brown     test: # par load, par interpolation, parmetis repartitioner
1770c4762a1bSJed Brown       suffix: 7_par_hdf5_int_parmetis
1771c4762a1bSJed Brown       requires: parmetis
1772c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1773c4762a1bSJed Brown       args: -interpolate after_create
1774c4762a1bSJed Brown     test: # par load, parmetis partitioner, par interpolation
1775c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1776c4762a1bSJed Brown       suffix: 7_par_hdf5_parmetis_int
1777c4762a1bSJed Brown       requires: parmetis
1778c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1779c4762a1bSJed Brown       args: -interpolate after_distribute
1780c4762a1bSJed Brown 
1781c4762a1bSJed Brown     test:
1782c4762a1bSJed Brown       suffix: 7_hdf5_hierarch
1783c4762a1bSJed Brown       requires: hdf5 ptscotch !complex
1784c4762a1bSJed Brown       nsize: {{2 3 4}separate output}
1785c4762a1bSJed Brown       args: -distribute
1786c4762a1bSJed Brown       args: -interpolate after_create
1787c4762a1bSJed Brown       args: -petscpartitioner_type matpartitioning -petscpartitioner_view ::ascii_info
1788c4762a1bSJed Brown       args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 2
1789c4762a1bSJed Brown       args: -mat_partitioning_hierarchical_coarseparttype ptscotch -mat_partitioning_hierarchical_fineparttype ptscotch
1790c4762a1bSJed Brown 
1791c4762a1bSJed Brown   test:
1792c4762a1bSJed Brown     suffix: 8
1793c4762a1bSJed Brown     requires: hdf5 !complex
1794c4762a1bSJed Brown     nsize: 4
1795c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1796c4762a1bSJed Brown     args: -distribute 0 -interpolate after_create
1797c4762a1bSJed Brown     args: -view_vertices_from_coords 0.,1.,0.,-0.5,1.,0.,0.583,-0.644,0.,-2.,-2.,-2. -view_vertices_from_coords_tol 1e-3
1798c4762a1bSJed Brown     args: -dm_plex_check_all
1799c4762a1bSJed Brown     args: -custom_view
1800c4762a1bSJed Brown 
1801c4762a1bSJed Brown   testset: # 9 HDF5 SEQUANTIAL LOAD
1802c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1803c4762a1bSJed Brown     args: -dm_plex_check_all
1804c4762a1bSJed Brown     args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates
1805c4762a1bSJed Brown     args: -dm_plex_hdf5_force_sequential
1806c4762a1bSJed Brown     args: -distribute
1807c4762a1bSJed Brown     test: # seq load, simple partitioner
1808c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple
1809c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1810c4762a1bSJed Brown       args: -interpolate none
1811c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1812c4762a1bSJed Brown       suffix: 9_seq_hdf5_int_simple
1813c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1814c4762a1bSJed Brown       args: -interpolate after_create
1815c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1816c4762a1bSJed Brown       nsize: {{2 4 5}}
1817c4762a1bSJed Brown       suffix: 9_seq_hdf5_int_metis
1818c4762a1bSJed Brown       requires: parmetis
1819c4762a1bSJed Brown       args: -interpolate after_create
1820c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1821c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1822c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple_int
1823c4762a1bSJed Brown       nsize: {{2 4 5}}
1824c4762a1bSJed Brown       args: -interpolate after_distribute
1825c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1826c4762a1bSJed Brown       # This is like 9_seq_hdf5_simple_int but testing error output of DMPlexCheckPointSFHeavy().
1827c4762a1bSJed Brown       # Once 9_seq_hdf5_simple_int gets fixed, this one gets broken.
1828c4762a1bSJed Brown       # We can then provide an intentionally broken mesh instead.
1829c4762a1bSJed Brown       TODO: This test is broken because PointSF is fixed.
1830c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple_int_err
1831c4762a1bSJed Brown       nsize: 4
1832c4762a1bSJed Brown       args: -interpolate after_distribute
1833c4762a1bSJed Brown       filter: sed -e "/PETSC ERROR/,$$d"
1834c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1835c4762a1bSJed Brown       nsize: {{2 4 5}}
1836c4762a1bSJed Brown       suffix: 9_seq_hdf5_metis_int
1837c4762a1bSJed Brown       requires: parmetis
1838c4762a1bSJed Brown       args: -interpolate after_distribute
1839c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1840c4762a1bSJed Brown 
1841c4762a1bSJed Brown   testset: # 9 HDF5 PARALLEL LOAD
1842c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1843c4762a1bSJed Brown     nsize: {{2 4 5}}
1844c4762a1bSJed Brown     args: -dm_plex_check_all
1845c4762a1bSJed Brown     args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /cells -dm_plex_hdf5_geometry_path /coordinates
1846c4762a1bSJed Brown     test: # par load
1847c4762a1bSJed Brown       suffix: 9_par_hdf5
1848c4762a1bSJed Brown       args: -interpolate none
1849c4762a1bSJed Brown     test: # par load, par interpolation
1850c4762a1bSJed Brown       suffix: 9_par_hdf5_int
1851c4762a1bSJed Brown       args: -interpolate after_create
1852c4762a1bSJed Brown     test: # par load, parmetis repartitioner
1853c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1854c4762a1bSJed Brown       suffix: 9_par_hdf5_parmetis
1855c4762a1bSJed Brown       requires: parmetis
1856c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1857c4762a1bSJed Brown       args: -interpolate none
1858c4762a1bSJed Brown     test: # par load, par interpolation, parmetis repartitioner
1859c4762a1bSJed Brown       suffix: 9_par_hdf5_int_parmetis
1860c4762a1bSJed Brown       requires: parmetis
1861c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1862c4762a1bSJed Brown       args: -interpolate after_create
1863c4762a1bSJed Brown     test: # par load, parmetis partitioner, par interpolation
1864c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1865c4762a1bSJed Brown       suffix: 9_par_hdf5_parmetis_int
1866c4762a1bSJed Brown       requires: parmetis
1867c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1868c4762a1bSJed Brown       args: -interpolate after_distribute
1869c4762a1bSJed Brown 
1870c4762a1bSJed Brown   testset: # 10 HDF5 PARALLEL LOAD
1871c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1872c4762a1bSJed Brown     nsize: {{2 4 7}}
1873c4762a1bSJed Brown     args: -dm_plex_check_all
1874c4762a1bSJed Brown     args: -filename ${DATAFILESPATH}/meshes/cube-hexahedra-refined2.h5 -dm_plex_create_from_hdf5_xdmf -dm_plex_hdf5_topology_path /topo -dm_plex_hdf5_geometry_path /geom
1875c4762a1bSJed Brown     test: # par load, par interpolation
1876c4762a1bSJed Brown       suffix: 10_par_hdf5_int
1877c4762a1bSJed Brown       args: -interpolate after_create
1878c4762a1bSJed Brown TEST*/
1879