xref: /petsc/src/dm/impls/plex/tests/ex18.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests for parallel mesh loading and parallel topological interpolation\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petsc/private/dmpleximpl.h>
4*c4762a1bSJed Brown /* List of test meshes
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown Network
7*c4762a1bSJed Brown -------
8*c4762a1bSJed Brown Test 0 (2 ranks):
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown network=0:
11*c4762a1bSJed Brown ---------
12*c4762a1bSJed Brown   cell 0   cell 1   cell 2          nCells-1       (edge)
13*c4762a1bSJed Brown 0 ------ 1 ------ 2 ------ 3 -- -- v --  -- nCells (vertex)
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   vertex distribution:
16*c4762a1bSJed Brown     rank 0: 0 1
17*c4762a1bSJed Brown     rank 1: 2 3 ... nCells
18*c4762a1bSJed Brown   cell(edge) distribution:
19*c4762a1bSJed Brown     rank 0: 0 1
20*c4762a1bSJed Brown     rank 1: 2 ... nCells-1
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown network=1:
23*c4762a1bSJed Brown ---------
24*c4762a1bSJed Brown                v2
25*c4762a1bSJed Brown                 ^
26*c4762a1bSJed Brown                 |
27*c4762a1bSJed Brown                cell 2
28*c4762a1bSJed Brown                 |
29*c4762a1bSJed Brown  v0 --cell 0--> v3--cell 1--> v1
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown   vertex distribution:
32*c4762a1bSJed Brown     rank 0: 0 1 3
33*c4762a1bSJed Brown     rank 1: 2
34*c4762a1bSJed Brown   cell(edge) distribution:
35*c4762a1bSJed Brown     rank 0: 0 1
36*c4762a1bSJed Brown     rank 1: 2
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   example:
39*c4762a1bSJed Brown     mpiexec -n 2 ./ex18 -distribute 1 -dim 1 -orig_dm_view -dist_dm_view -dist_dm_view -petscpartitioner_type parmetis -ncells 50
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown Triangle
42*c4762a1bSJed Brown --------
43*c4762a1bSJed Brown Test 0 (2 ranks):
44*c4762a1bSJed Brown Two triangles sharing a face
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown         2
47*c4762a1bSJed Brown       / | \
48*c4762a1bSJed Brown      /  |  \
49*c4762a1bSJed Brown     /   |   \
50*c4762a1bSJed Brown    0  0 | 1  3
51*c4762a1bSJed Brown     \   |   /
52*c4762a1bSJed Brown      \  |  /
53*c4762a1bSJed Brown       \ | /
54*c4762a1bSJed Brown         1
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown   vertex distribution:
57*c4762a1bSJed Brown     rank 0: 0 1
58*c4762a1bSJed Brown     rank 1: 2 3
59*c4762a1bSJed Brown   cell distribution:
60*c4762a1bSJed Brown     rank 0: 0
61*c4762a1bSJed Brown     rank 1: 1
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown Test 1 (3 ranks):
64*c4762a1bSJed Brown Four triangles partitioned across 3 ranks
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown    0 _______ 3
67*c4762a1bSJed Brown    | \     / |
68*c4762a1bSJed Brown    |  \ 1 /  |
69*c4762a1bSJed Brown    |   \ /   |
70*c4762a1bSJed Brown    | 0  2  2 |
71*c4762a1bSJed Brown    |   / \   |
72*c4762a1bSJed Brown    |  / 3 \  |
73*c4762a1bSJed Brown    | /     \ |
74*c4762a1bSJed Brown    1 ------- 4
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   vertex distribution:
77*c4762a1bSJed Brown     rank 0: 0 1
78*c4762a1bSJed Brown     rank 1: 2 3
79*c4762a1bSJed Brown     rank 2: 4
80*c4762a1bSJed Brown   cell distribution:
81*c4762a1bSJed Brown     rank 0: 0
82*c4762a1bSJed Brown     rank 1: 1
83*c4762a1bSJed Brown     rank 2: 2 3
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown Test 2 (3 ranks):
86*c4762a1bSJed Brown Four triangles partitioned across 3 ranks
87*c4762a1bSJed Brown 
88*c4762a1bSJed Brown    1 _______ 3
89*c4762a1bSJed Brown    | \     / |
90*c4762a1bSJed Brown    |  \ 1 /  |
91*c4762a1bSJed Brown    |   \ /   |
92*c4762a1bSJed Brown    | 0  0  2 |
93*c4762a1bSJed Brown    |   / \   |
94*c4762a1bSJed Brown    |  / 3 \  |
95*c4762a1bSJed Brown    | /     \ |
96*c4762a1bSJed Brown    2 ------- 4
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   vertex distribution:
99*c4762a1bSJed Brown     rank 0: 0 1
100*c4762a1bSJed Brown     rank 1: 2 3
101*c4762a1bSJed Brown     rank 2: 4
102*c4762a1bSJed Brown   cell distribution:
103*c4762a1bSJed Brown     rank 0: 0
104*c4762a1bSJed Brown     rank 1: 1
105*c4762a1bSJed Brown     rank 2: 2 3
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown Tetrahedron
108*c4762a1bSJed Brown -----------
109*c4762a1bSJed Brown Test 0:
110*c4762a1bSJed Brown Two tets sharing a face
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown  cell   3 _______    cell
113*c4762a1bSJed Brown  0    / | \      \   1
114*c4762a1bSJed Brown      /  |  \      \
115*c4762a1bSJed Brown     /   |   \      \
116*c4762a1bSJed Brown    0----|----4-----2
117*c4762a1bSJed Brown     \   |   /      /
118*c4762a1bSJed Brown      \  |  /      /
119*c4762a1bSJed Brown       \ | /      /
120*c4762a1bSJed Brown         1-------
121*c4762a1bSJed Brown    y
122*c4762a1bSJed Brown    | x
123*c4762a1bSJed Brown    |/
124*c4762a1bSJed Brown    *----z
125*c4762a1bSJed Brown 
126*c4762a1bSJed Brown   vertex distribution:
127*c4762a1bSJed Brown     rank 0: 0 1
128*c4762a1bSJed Brown     rank 1: 2 3 4
129*c4762a1bSJed Brown   cell distribution:
130*c4762a1bSJed Brown     rank 0: 0
131*c4762a1bSJed Brown     rank 1: 1
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown Quadrilateral
134*c4762a1bSJed Brown -------------
135*c4762a1bSJed Brown Test 0 (2 ranks):
136*c4762a1bSJed Brown Two quads sharing a face
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown    3-------2-------5
139*c4762a1bSJed Brown    |       |       |
140*c4762a1bSJed Brown    |   0   |   1   |
141*c4762a1bSJed Brown    |       |       |
142*c4762a1bSJed Brown    0-------1-------4
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown   vertex distribution:
145*c4762a1bSJed Brown     rank 0: 0 1 2
146*c4762a1bSJed Brown     rank 1: 3 4 5
147*c4762a1bSJed Brown   cell distribution:
148*c4762a1bSJed Brown     rank 0: 0
149*c4762a1bSJed Brown     rank 1: 1
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown TODO Test 1:
152*c4762a1bSJed Brown A quad and a triangle sharing a face
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown    5-------4
155*c4762a1bSJed Brown    |       | \
156*c4762a1bSJed Brown    |   0   |  \
157*c4762a1bSJed Brown    |       | 1 \
158*c4762a1bSJed Brown    2-------3----6
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown Hexahedron
161*c4762a1bSJed Brown ----------
162*c4762a1bSJed Brown Test 0 (2 ranks):
163*c4762a1bSJed Brown Two hexes sharing a face
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown cell   7-------------6-------------11 cell
166*c4762a1bSJed Brown 0     /|            /|            /|     1
167*c4762a1bSJed Brown      / |   F1      / |   F7      / |
168*c4762a1bSJed Brown     /  |          /  |          /  |
169*c4762a1bSJed Brown    4-------------5-------------10  |
170*c4762a1bSJed Brown    |   |     F4  |   |     F10 |   |
171*c4762a1bSJed Brown    |   |         |   |         |   |
172*c4762a1bSJed Brown    |F5 |         |F3 |         |F9 |
173*c4762a1bSJed Brown    |   |  F2     |   |   F8    |   |
174*c4762a1bSJed Brown    |   3---------|---2---------|---9
175*c4762a1bSJed Brown    |  /          |  /          |  /
176*c4762a1bSJed Brown    | /   F0      | /    F6     | /
177*c4762a1bSJed Brown    |/            |/            |/
178*c4762a1bSJed Brown    0-------------1-------------8
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   vertex distribution:
181*c4762a1bSJed Brown     rank 0: 0 1 2 3 4 5
182*c4762a1bSJed Brown     rank 1: 6 7 8 9 10 11
183*c4762a1bSJed Brown   cell distribution:
184*c4762a1bSJed Brown     rank 0: 0
185*c4762a1bSJed Brown     rank 1: 1
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown */
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown typedef enum {NONE, CREATE, AFTER_CREATE, AFTER_DISTRIBUTE} InterpType;
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown typedef struct {
192*c4762a1bSJed Brown   PetscInt   debug;                        /* The debugging level */
193*c4762a1bSJed Brown   PetscInt   testNum;                      /* Indicates the mesh to create */
194*c4762a1bSJed Brown   PetscInt   dim;                          /* The topological mesh dimension */
195*c4762a1bSJed Brown   PetscBool  cellSimplex;                  /* Use simplices or hexes */
196*c4762a1bSJed Brown   PetscBool  distribute;                   /* Distribute the mesh */
197*c4762a1bSJed Brown   InterpType interpolate;                  /* Interpolate the mesh before or after DMPlexDistribute() */
198*c4762a1bSJed Brown   PetscBool  useGenerator;                 /* Construct mesh with a mesh generator */
199*c4762a1bSJed Brown   PetscBool  testOrientIF;                 /* Test for different original interface orientations */
200*c4762a1bSJed Brown   PetscBool  testHeavy;                    /* Run the heavy PointSF test */
201*c4762a1bSJed Brown   PetscBool  customView;                   /* Show results of DMPlexIsInterpolated() etc. */
202*c4762a1bSJed Brown   PetscInt   ornt[2];                      /* Orientation of interface on rank 0 and rank 1 */
203*c4762a1bSJed Brown   PetscInt   faces[3];                     /* Number of faces per dimension for generator */
204*c4762a1bSJed Brown   PetscReal  coords[128];
205*c4762a1bSJed Brown   PetscReal  coordsTol;
206*c4762a1bSJed Brown   PetscInt   ncoords;
207*c4762a1bSJed Brown   PetscInt   pointsToExpand[128];
208*c4762a1bSJed Brown   PetscInt   nPointsToExpand;
209*c4762a1bSJed Brown   PetscBool  testExpandPointsEmpty;
210*c4762a1bSJed Brown   char       filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
211*c4762a1bSJed Brown } AppCtx;
212*c4762a1bSJed Brown 
213*c4762a1bSJed Brown struct _n_PortableBoundary {
214*c4762a1bSJed Brown   Vec coordinates;
215*c4762a1bSJed Brown   PetscInt depth;
216*c4762a1bSJed Brown   PetscSection *sections;
217*c4762a1bSJed Brown };
218*c4762a1bSJed Brown typedef struct _n_PortableBoundary * PortableBoundary;
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown static PetscLogStage  stage[3];
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM, PortableBoundary);
223*c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM,PetscBool);
224*c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM, PortableBoundary *);
225*c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM, IS, PetscSection, IS *);
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown static PetscErrorCode PortableBoundaryDestroy(PortableBoundary *bnd)
228*c4762a1bSJed Brown {
229*c4762a1bSJed Brown   PetscInt       d;
230*c4762a1bSJed Brown   PetscErrorCode ierr;
231*c4762a1bSJed Brown 
232*c4762a1bSJed Brown   PetscFunctionBegin;
233*c4762a1bSJed Brown   if (!*bnd) PetscFunctionReturn(0);
234*c4762a1bSJed Brown   ierr = VecDestroy(&(*bnd)->coordinates);CHKERRQ(ierr);
235*c4762a1bSJed Brown   for (d=0; d < (*bnd)->depth; d++) {
236*c4762a1bSJed Brown     ierr = PetscSectionDestroy(&(*bnd)->sections[d]);CHKERRQ(ierr);
237*c4762a1bSJed Brown   }
238*c4762a1bSJed Brown   ierr = PetscFree((*bnd)->sections);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = PetscFree(*bnd);CHKERRQ(ierr);
240*c4762a1bSJed Brown   PetscFunctionReturn(0);
241*c4762a1bSJed Brown }
242*c4762a1bSJed Brown 
243*c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
244*c4762a1bSJed Brown {
245*c4762a1bSJed Brown   const char    *interpTypes[4]  = {"none", "create", "after_create", "after_distribute"};
246*c4762a1bSJed Brown   PetscInt       interp=NONE, dim;
247*c4762a1bSJed Brown   PetscBool      flg1, flg2;
248*c4762a1bSJed Brown   PetscErrorCode ierr;
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown   PetscFunctionBegin;
251*c4762a1bSJed Brown   options->debug        = 0;
252*c4762a1bSJed Brown   options->testNum      = 0;
253*c4762a1bSJed Brown   options->dim          = 2;
254*c4762a1bSJed Brown   options->cellSimplex  = PETSC_TRUE;
255*c4762a1bSJed Brown   options->distribute   = PETSC_FALSE;
256*c4762a1bSJed Brown   options->interpolate  = NONE;
257*c4762a1bSJed Brown   options->useGenerator = PETSC_FALSE;
258*c4762a1bSJed Brown   options->testOrientIF = PETSC_FALSE;
259*c4762a1bSJed Brown   options->testHeavy    = PETSC_TRUE;
260*c4762a1bSJed Brown   options->customView   = PETSC_FALSE;
261*c4762a1bSJed Brown   options->testExpandPointsEmpty = PETSC_FALSE;
262*c4762a1bSJed Brown   options->ornt[0]      = 0;
263*c4762a1bSJed Brown   options->ornt[1]      = 0;
264*c4762a1bSJed Brown   options->faces[0]     = 2;
265*c4762a1bSJed Brown   options->faces[1]     = 2;
266*c4762a1bSJed Brown   options->faces[2]     = 2;
267*c4762a1bSJed Brown   options->filename[0]  = '\0';
268*c4762a1bSJed Brown   options->coordsTol    = PETSC_DEFAULT;
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex18.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex18.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Generate simplices if true, otherwise hexes", "ex18.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = PetscOptionsBool("-distribute", "Distribute the mesh", "ex18.c", options->distribute, &options->distribute, NULL);CHKERRQ(ierr);
275*c4762a1bSJed Brown   ierr = PetscOptionsEList("-interpolate", "Type of mesh interpolation (none, create, after_create, after_distribute)", "ex18.c", interpTypes, 4, interpTypes[options->interpolate], &interp, NULL);CHKERRQ(ierr);
276*c4762a1bSJed Brown   options->interpolate = (InterpType) interp;
277*c4762a1bSJed Brown   if (!options->distribute && options->interpolate == AFTER_DISTRIBUTE) SETERRQ(comm, PETSC_ERR_SUP, "-interpolate after_distribute  needs  -distribute 1");
278*c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex18.c", options->useGenerator, &options->useGenerator, NULL);CHKERRQ(ierr);
279*c4762a1bSJed Brown   options->ncoords = 128;
280*c4762a1bSJed Brown   ierr = PetscOptionsRealArray("-view_vertices_from_coords", "Print DAG points corresponding to vertices with given coordinates", "ex18.c", options->coords, &options->ncoords, NULL);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = PetscOptionsReal("-view_vertices_from_coords_tol", "Tolerance for -view_vertices_from_coords", "ex18.c", options->coordsTol, &options->coordsTol, NULL);CHKERRQ(ierr);
282*c4762a1bSJed Brown   options->nPointsToExpand = 128;
283*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-test_expand_points", "Expand given array of DAG point using DMPlexGetConeRecursive() and print results", "ex18.c", options->pointsToExpand, &options->nPointsToExpand, NULL);CHKERRQ(ierr);
284*c4762a1bSJed Brown   if (options->nPointsToExpand) {
285*c4762a1bSJed Brown     ierr = PetscOptionsBool("-test_expand_points_empty", "For -test_expand_points, rank 0 will have empty input array", "ex18.c", options->testExpandPointsEmpty, &options->testExpandPointsEmpty, NULL);CHKERRQ(ierr);
286*c4762a1bSJed Brown   }
287*c4762a1bSJed Brown   ierr = PetscOptionsBool("-test_heavy", "Run the heavy PointSF test", "ex18.c", options->testHeavy, &options->testHeavy, NULL);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = PetscOptionsBool("-custom_view", "Custom DMPlex view", "ex18.c", options->customView, &options->customView, NULL);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex18.c", options->dim, &options->dim, &flg1,1,3);CHKERRQ(ierr);
290*c4762a1bSJed Brown   dim = 3;
291*c4762a1bSJed Brown   ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex18.c", options->faces, &dim, &flg2);CHKERRQ(ierr);
292*c4762a1bSJed Brown   if (flg2) {
293*c4762a1bSJed Brown     if (flg1 && dim != options->dim) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "specified -dim %D is not equal to length %D of -faces (note that -dim can be omitted)", options->dim, dim);
294*c4762a1bSJed Brown     options->dim = dim;
295*c4762a1bSJed Brown   }
296*c4762a1bSJed Brown   ierr = PetscOptionsString("-filename", "The mesh file", "ex18.c", options->filename, options->filename, PETSC_MAX_PATH_LEN, NULL);CHKERRQ(ierr);
297*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
298*c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
299*c4762a1bSJed Brown   if (flg2 != options->testOrientIF) SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "neither or both -rotate_interface_0 -rotate_interface_1 must be set");
300*c4762a1bSJed Brown   if (options->testOrientIF) {
301*c4762a1bSJed Brown     PetscInt i;
302*c4762a1bSJed Brown     for (i=0; i<2; i++) {
303*c4762a1bSJed Brown       if (options->ornt[i] >= 10) options->ornt[i] = -(options->ornt[i]-10);  /* 11 12 13 become -1 -2 -3 */
304*c4762a1bSJed Brown     }
305*c4762a1bSJed Brown     options->filename[0]  = 0;
306*c4762a1bSJed Brown     options->useGenerator = PETSC_FALSE;
307*c4762a1bSJed Brown     options->dim          = 3;
308*c4762a1bSJed Brown     options->cellSimplex  = PETSC_TRUE;
309*c4762a1bSJed Brown     options->interpolate  = CREATE;
310*c4762a1bSJed Brown     options->distribute   = PETSC_FALSE;
311*c4762a1bSJed Brown   }
312*c4762a1bSJed Brown   ierr = PetscOptionsEnd();
313*c4762a1bSJed Brown   PetscFunctionReturn(0);
314*c4762a1bSJed Brown }
315*c4762a1bSJed Brown 
316*c4762a1bSJed Brown static PetscErrorCode CreateMesh_1D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
317*c4762a1bSJed Brown {
318*c4762a1bSJed Brown   PetscInt       testNum = user->testNum;
319*c4762a1bSJed Brown   PetscMPIInt    rank,size;
320*c4762a1bSJed Brown   PetscErrorCode ierr;
321*c4762a1bSJed Brown   PetscInt       numCorners=2,i;
322*c4762a1bSJed Brown   PetscInt       numCells,numVertices,network;
323*c4762a1bSJed Brown   int            *cells;
324*c4762a1bSJed Brown   PetscReal      *coords;
325*c4762a1bSJed Brown 
326*c4762a1bSJed Brown   PetscFunctionBegin;
327*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
328*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
329*c4762a1bSJed Brown   if (size > 2) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for <=2 processes",testNum);
330*c4762a1bSJed Brown 
331*c4762a1bSJed Brown   numCells = 3;
332*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL, NULL, "-ncells", &numCells, NULL);CHKERRQ(ierr);
333*c4762a1bSJed Brown   if (numCells < 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test ncells must >=3",numCells);
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown   if (size == 1) {
336*c4762a1bSJed Brown     double *dcoords;
337*c4762a1bSJed Brown     numVertices = numCells + 1;
338*c4762a1bSJed Brown     ierr = PetscMalloc2(2*numCells,&cells,2*numVertices,&dcoords);CHKERRQ(ierr);
339*c4762a1bSJed Brown     for (i=0; i<numCells; i++) {
340*c4762a1bSJed Brown       cells[2*i] = i; cells[2*i+1] = i + 1;
341*c4762a1bSJed Brown       dcoords[2*i] = i; dcoords[2*i+1] = i + 1;
342*c4762a1bSJed Brown     }
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown     ierr = DMPlexCreateFromCellList(comm, user->dim, numCells, numVertices, numCorners, PETSC_FALSE, cells, user->dim, dcoords, dm);CHKERRQ(ierr);
345*c4762a1bSJed Brown     ierr = PetscFree2(cells,coords);CHKERRQ(ierr);
346*c4762a1bSJed Brown     PetscFunctionReturn(0);
347*c4762a1bSJed Brown   }
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown   network = 0;
350*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL, NULL, "-network_case", &network, NULL);CHKERRQ(ierr);
351*c4762a1bSJed Brown   if (network == 0) {
352*c4762a1bSJed Brown     switch (rank) {
353*c4762a1bSJed Brown     case 0:
354*c4762a1bSJed Brown     {
355*c4762a1bSJed Brown       numCells    = 2;
356*c4762a1bSJed Brown       numVertices = numCells;
357*c4762a1bSJed Brown       ierr = PetscMalloc2(2*numCells,&cells,2*numCells,&coords);CHKERRQ(ierr);
358*c4762a1bSJed Brown       cells[0] = 0; cells[1] = 1;
359*c4762a1bSJed Brown       cells[2] = 1; cells[3] = 2;
360*c4762a1bSJed Brown       coords[0] = 0.; coords[1] = 1.;
361*c4762a1bSJed Brown       coords[2] = 1.; coords[3] = 2.;
362*c4762a1bSJed Brown     }
363*c4762a1bSJed Brown     break;
364*c4762a1bSJed Brown     case 1:
365*c4762a1bSJed Brown     {
366*c4762a1bSJed Brown       numCells    -= 2;
367*c4762a1bSJed Brown       numVertices = numCells + 1;
368*c4762a1bSJed Brown       ierr = PetscMalloc2(2*numCells,&cells,2*numCells,&coords);CHKERRQ(ierr);
369*c4762a1bSJed Brown       for (i=0; i<numCells; i++) {
370*c4762a1bSJed Brown         cells[2*i] = 2+i; cells[2*i+1] = 2 + i + 1;
371*c4762a1bSJed Brown         coords[2*i] = 2+i; coords[2*i+1] = 2 + i + 1;
372*c4762a1bSJed Brown       }
373*c4762a1bSJed Brown     }
374*c4762a1bSJed Brown     break;
375*c4762a1bSJed Brown     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
376*c4762a1bSJed Brown     }
377*c4762a1bSJed Brown   } else { /* network_case = 1 */
378*c4762a1bSJed Brown     /* ----------------------- */
379*c4762a1bSJed Brown     switch (rank) {
380*c4762a1bSJed Brown     case 0:
381*c4762a1bSJed Brown     {
382*c4762a1bSJed Brown       numCells    = 2;
383*c4762a1bSJed Brown       numVertices = 3;
384*c4762a1bSJed Brown       ierr = PetscMalloc2(2*numCells,&cells,2*numCells,&coords);CHKERRQ(ierr);
385*c4762a1bSJed Brown       cells[0] = 0; cells[1] = 3;
386*c4762a1bSJed Brown       cells[2] = 3; cells[3] = 1;
387*c4762a1bSJed Brown     }
388*c4762a1bSJed Brown     break;
389*c4762a1bSJed Brown     case 1:
390*c4762a1bSJed Brown     {
391*c4762a1bSJed Brown       numCells    = 1;
392*c4762a1bSJed Brown       numVertices = 1;
393*c4762a1bSJed Brown       ierr = PetscMalloc2(2*numCells,&cells,2*numCells,&coords);CHKERRQ(ierr);
394*c4762a1bSJed Brown       cells[0] = 3; cells[1] = 2;
395*c4762a1bSJed Brown     }
396*c4762a1bSJed Brown     break;
397*c4762a1bSJed Brown     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
398*c4762a1bSJed Brown     }
399*c4762a1bSJed Brown   }
400*c4762a1bSJed Brown   ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, PETSC_FALSE, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
401*c4762a1bSJed Brown   ierr = PetscFree2(cells,coords);CHKERRQ(ierr);
402*c4762a1bSJed Brown   PetscFunctionReturn(0);
403*c4762a1bSJed Brown }
404*c4762a1bSJed Brown 
405*c4762a1bSJed Brown static PetscErrorCode CreateSimplex_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
406*c4762a1bSJed Brown {
407*c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
408*c4762a1bSJed Brown   PetscMPIInt    rank, size;
409*c4762a1bSJed Brown   PetscErrorCode ierr;
410*c4762a1bSJed Brown 
411*c4762a1bSJed Brown   PetscFunctionBegin;
412*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
413*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
414*c4762a1bSJed Brown   switch (testNum) {
415*c4762a1bSJed Brown   case 0:
416*c4762a1bSJed Brown     if (size != 2) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 2 processes", testNum);
417*c4762a1bSJed Brown     switch (rank) {
418*c4762a1bSJed Brown       case 0:
419*c4762a1bSJed Brown       {
420*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
421*c4762a1bSJed Brown         const int      cells[3]  = {0, 1, 2};
422*c4762a1bSJed Brown         PetscReal      coords[4] = {-0.5, 0.5, 0.0, 0.0};
423*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
424*c4762a1bSJed Brown 
425*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
426*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
427*c4762a1bSJed Brown       }
428*c4762a1bSJed Brown       break;
429*c4762a1bSJed Brown       case 1:
430*c4762a1bSJed Brown       {
431*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
432*c4762a1bSJed Brown         const int      cells[3]  = {1, 3, 2};
433*c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 1.0, 0.5, 0.5};
434*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
435*c4762a1bSJed Brown 
436*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
437*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
438*c4762a1bSJed Brown       }
439*c4762a1bSJed Brown       break;
440*c4762a1bSJed Brown       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
441*c4762a1bSJed Brown     }
442*c4762a1bSJed Brown     break;
443*c4762a1bSJed Brown   case 1:
444*c4762a1bSJed Brown     if (size != 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 3 processes", testNum);
445*c4762a1bSJed Brown     switch (rank) {
446*c4762a1bSJed Brown       case 0:
447*c4762a1bSJed Brown       {
448*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
449*c4762a1bSJed Brown         const int      cells[3]  = {0, 1, 2};
450*c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 1.0, 0.0, 0.0};
451*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
452*c4762a1bSJed Brown 
453*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
454*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
455*c4762a1bSJed Brown       }
456*c4762a1bSJed Brown       break;
457*c4762a1bSJed Brown       case 1:
458*c4762a1bSJed Brown       {
459*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
460*c4762a1bSJed Brown         const int      cells[3]  = {0, 2, 3};
461*c4762a1bSJed Brown         PetscReal      coords[4] = {0.5, 0.5, 1.0, 1.0};
462*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
463*c4762a1bSJed Brown 
464*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
465*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
466*c4762a1bSJed Brown       }
467*c4762a1bSJed Brown       break;
468*c4762a1bSJed Brown       case 2:
469*c4762a1bSJed Brown       {
470*c4762a1bSJed Brown         const PetscInt numCells  = 2, numVertices = 1, numCorners = 3;
471*c4762a1bSJed Brown         const int      cells[6]  = {2, 4, 3, 2, 1, 4};
472*c4762a1bSJed Brown         PetscReal      coords[2] = {1.0, 0.0};
473*c4762a1bSJed Brown         PetscInt       markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
474*c4762a1bSJed Brown 
475*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
476*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
477*c4762a1bSJed Brown       }
478*c4762a1bSJed Brown       break;
479*c4762a1bSJed Brown       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
480*c4762a1bSJed Brown     }
481*c4762a1bSJed Brown     break;
482*c4762a1bSJed Brown   case 2:
483*c4762a1bSJed Brown     if (size != 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 3 processes", testNum);
484*c4762a1bSJed Brown     switch (rank) {
485*c4762a1bSJed Brown       case 0:
486*c4762a1bSJed Brown       {
487*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
488*c4762a1bSJed Brown         const int      cells[3]  = {1, 2, 0};
489*c4762a1bSJed Brown         PetscReal      coords[4] = {0.5, 0.5, 0.0, 1.0};
490*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
491*c4762a1bSJed Brown 
492*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
493*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
494*c4762a1bSJed Brown       }
495*c4762a1bSJed Brown       break;
496*c4762a1bSJed Brown       case 1:
497*c4762a1bSJed Brown       {
498*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 3;
499*c4762a1bSJed Brown         const int      cells[3]  = {1, 0, 3};
500*c4762a1bSJed Brown         PetscReal      coords[4] = {0.0, 0.0, 1.0, 1.0};
501*c4762a1bSJed Brown         PetscInt       markerPoints[6] = {1, 1, 2, 1, 3, 1};
502*c4762a1bSJed Brown 
503*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
504*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
505*c4762a1bSJed Brown       }
506*c4762a1bSJed Brown       break;
507*c4762a1bSJed Brown       case 2:
508*c4762a1bSJed Brown       {
509*c4762a1bSJed Brown         const PetscInt numCells  = 2, numVertices = 1, numCorners = 3;
510*c4762a1bSJed Brown         const int      cells[6]  = {0, 4, 3, 0, 2, 4};
511*c4762a1bSJed Brown         PetscReal      coords[2] = {1.0, 0.0};
512*c4762a1bSJed Brown         PetscInt       markerPoints[10] = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
513*c4762a1bSJed Brown 
514*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
515*c4762a1bSJed Brown         for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
516*c4762a1bSJed Brown       }
517*c4762a1bSJed Brown       break;
518*c4762a1bSJed Brown       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
519*c4762a1bSJed Brown     }
520*c4762a1bSJed Brown     break;
521*c4762a1bSJed Brown   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %D", testNum);
522*c4762a1bSJed Brown   }
523*c4762a1bSJed Brown   PetscFunctionReturn(0);
524*c4762a1bSJed Brown }
525*c4762a1bSJed Brown 
526*c4762a1bSJed Brown static PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
527*c4762a1bSJed Brown {
528*c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
529*c4762a1bSJed Brown   PetscMPIInt    rank, size;
530*c4762a1bSJed Brown   PetscErrorCode ierr;
531*c4762a1bSJed Brown 
532*c4762a1bSJed Brown   PetscFunctionBegin;
533*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
534*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
535*c4762a1bSJed Brown   switch (testNum) {
536*c4762a1bSJed Brown   case 0:
537*c4762a1bSJed Brown     if (size != 2) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 2 processes", testNum);
538*c4762a1bSJed Brown     switch (rank) {
539*c4762a1bSJed Brown       case 0:
540*c4762a1bSJed Brown       {
541*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 2, numCorners = 4;
542*c4762a1bSJed Brown         const int      cells[4]  = {0, 2, 1, 3};
543*c4762a1bSJed Brown         PetscReal      coords[6] = {0.0, 0.0, -0.5,  0.0, -0.5, 0.0};
544*c4762a1bSJed Brown         PetscInt       markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1};
545*c4762a1bSJed Brown 
546*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
547*c4762a1bSJed Brown         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
548*c4762a1bSJed Brown       }
549*c4762a1bSJed Brown       break;
550*c4762a1bSJed Brown       case 1:
551*c4762a1bSJed Brown       {
552*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
553*c4762a1bSJed Brown         const int      cells[4]  = {1, 2, 4, 3};
554*c4762a1bSJed Brown         PetscReal      coords[9] = {1.0, 0.0, 0.0,  0.0, 0.5, 0.0,  0.0, 0.0, 0.5};
555*c4762a1bSJed Brown         PetscInt       markerPoints[8] = {1, 1, 2, 1, 3, 1, 4, 1};
556*c4762a1bSJed Brown 
557*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
558*c4762a1bSJed Brown         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
559*c4762a1bSJed Brown       }
560*c4762a1bSJed Brown       break;
561*c4762a1bSJed Brown       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
562*c4762a1bSJed Brown     }
563*c4762a1bSJed Brown     break;
564*c4762a1bSJed Brown   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %D", testNum);
565*c4762a1bSJed Brown   }
566*c4762a1bSJed Brown   if (user->testOrientIF) {
567*c4762a1bSJed Brown     PetscInt start;
568*c4762a1bSJed Brown     PetscBool reverse;
569*c4762a1bSJed Brown     PetscInt ifp[] = {8, 6};
570*c4762a1bSJed Brown 
571*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) *dm, "Mesh before orientation");CHKERRQ(ierr);
572*c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-before_orientation_dm_view");CHKERRQ(ierr);
573*c4762a1bSJed Brown     /* rotate interface face ifp[rank] by given orientation ornt[rank] */
574*c4762a1bSJed Brown     ierr = DMPlexFixFaceOrientations_Translate_Private(user->ornt[rank], &start, &reverse);CHKERRQ(ierr);
575*c4762a1bSJed Brown     ierr = DMPlexOrientCell_Internal(*dm, ifp[rank], start, reverse);CHKERRQ(ierr);
576*c4762a1bSJed Brown     ierr = DMPlexCheckFaces(*dm, 0);CHKERRQ(ierr);
577*c4762a1bSJed Brown     ierr = DMPlexOrientInterface_Internal(*dm);CHKERRQ(ierr);
578*c4762a1bSJed Brown     ierr = PetscPrintf(comm, "Orientation test PASSED\n");CHKERRQ(ierr);
579*c4762a1bSJed Brown   }
580*c4762a1bSJed Brown   PetscFunctionReturn(0);
581*c4762a1bSJed Brown }
582*c4762a1bSJed Brown 
583*c4762a1bSJed Brown static PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
584*c4762a1bSJed Brown {
585*c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
586*c4762a1bSJed Brown   PetscMPIInt    rank, size;
587*c4762a1bSJed Brown   PetscErrorCode ierr;
588*c4762a1bSJed Brown 
589*c4762a1bSJed Brown   PetscFunctionBegin;
590*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
591*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
592*c4762a1bSJed Brown   switch (testNum) {
593*c4762a1bSJed Brown   case 0:
594*c4762a1bSJed Brown     if (size != 2) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 2 processes", testNum);
595*c4762a1bSJed Brown     switch (rank) {
596*c4762a1bSJed Brown       case 0:
597*c4762a1bSJed Brown       {
598*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
599*c4762a1bSJed Brown         const int      cells[4]  = {0, 1, 2, 3};
600*c4762a1bSJed Brown         PetscReal      coords[6] = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0};
601*c4762a1bSJed Brown         PetscInt       markerPoints[4*2] = {1, 1, 2, 1, 3, 1, 4, 1};
602*c4762a1bSJed Brown 
603*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
604*c4762a1bSJed Brown         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
605*c4762a1bSJed Brown       }
606*c4762a1bSJed Brown       break;
607*c4762a1bSJed Brown       case 1:
608*c4762a1bSJed Brown       {
609*c4762a1bSJed Brown         const PetscInt numCells  = 1, numVertices = 3, numCorners = 4;
610*c4762a1bSJed Brown         const int      cells[4]  = {1, 4, 5, 2};
611*c4762a1bSJed Brown         PetscReal      coords[6] = {-0.5, 1.0, 0.5, 0.0, 0.5, 1.0};
612*c4762a1bSJed Brown         PetscInt       markerPoints[4*2] = {1, 1, 2, 1, 3, 1, 4, 1};
613*c4762a1bSJed Brown 
614*c4762a1bSJed Brown         ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
615*c4762a1bSJed Brown         for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
616*c4762a1bSJed Brown       }
617*c4762a1bSJed Brown       break;
618*c4762a1bSJed Brown       default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
619*c4762a1bSJed Brown     }
620*c4762a1bSJed Brown     break;
621*c4762a1bSJed Brown   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %D", testNum);
622*c4762a1bSJed Brown   }
623*c4762a1bSJed Brown   PetscFunctionReturn(0);
624*c4762a1bSJed Brown }
625*c4762a1bSJed Brown 
626*c4762a1bSJed Brown static PetscErrorCode CreateHex_3D(MPI_Comm comm, PetscBool interpolate, AppCtx *user, DM *dm)
627*c4762a1bSJed Brown {
628*c4762a1bSJed Brown   PetscInt       testNum = user->testNum, p;
629*c4762a1bSJed Brown   PetscMPIInt    rank, size;
630*c4762a1bSJed Brown   PetscErrorCode ierr;
631*c4762a1bSJed Brown 
632*c4762a1bSJed Brown   PetscFunctionBegin;
633*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
634*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
635*c4762a1bSJed Brown   switch (testNum) {
636*c4762a1bSJed Brown   case 0:
637*c4762a1bSJed Brown     if (size != 2) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Test mesh %d only for 2 processes", testNum);
638*c4762a1bSJed Brown     switch (rank) {
639*c4762a1bSJed Brown     case 0:
640*c4762a1bSJed Brown     {
641*c4762a1bSJed Brown       const PetscInt numCells  = 1, numVertices = 6, numCorners = 8;
642*c4762a1bSJed Brown       const int      cells[8]  = {0, 3, 2, 1, 4, 5, 6, 7};
643*c4762a1bSJed 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};
644*c4762a1bSJed Brown       PetscInt       markerPoints[8*2] = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
645*c4762a1bSJed Brown 
646*c4762a1bSJed Brown       ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
647*c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
648*c4762a1bSJed Brown     }
649*c4762a1bSJed Brown     break;
650*c4762a1bSJed Brown     case 1:
651*c4762a1bSJed Brown     {
652*c4762a1bSJed Brown       const PetscInt numCells  = 1, numVertices = 6, numCorners = 8;
653*c4762a1bSJed Brown       const int      cells[8]  = {1, 2, 9, 8, 5, 10, 11, 6};
654*c4762a1bSJed 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};
655*c4762a1bSJed Brown       PetscInt       markerPoints[8*2] = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
656*c4762a1bSJed Brown 
657*c4762a1bSJed Brown       ierr = DMPlexCreateFromCellListParallel(comm, user->dim, numCells, numVertices, numCorners, interpolate, cells, user->dim, coords, NULL, dm);CHKERRQ(ierr);
658*c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);}
659*c4762a1bSJed Brown     }
660*c4762a1bSJed Brown     break;
661*c4762a1bSJed Brown     default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh for rank %d", rank);
662*c4762a1bSJed Brown     }
663*c4762a1bSJed Brown   break;
664*c4762a1bSJed Brown   default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %D", testNum);
665*c4762a1bSJed Brown   }
666*c4762a1bSJed Brown   PetscFunctionReturn(0);
667*c4762a1bSJed Brown }
668*c4762a1bSJed Brown 
669*c4762a1bSJed Brown static PetscErrorCode CustomView(DM dm, PetscViewer v)
670*c4762a1bSJed Brown {
671*c4762a1bSJed Brown   DMPlexInterpolatedFlag interpolated;
672*c4762a1bSJed Brown   PetscBool              distributed;
673*c4762a1bSJed Brown   PetscErrorCode         ierr;
674*c4762a1bSJed Brown 
675*c4762a1bSJed Brown   PetscFunctionBegin;
676*c4762a1bSJed Brown   ierr = DMPlexIsDistributed(dm, &distributed);CHKERRQ(ierr);
677*c4762a1bSJed Brown   ierr = DMPlexIsInterpolatedCollective(dm, &interpolated);CHKERRQ(ierr);
678*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(v, "DMPlexIsDistributed: %s\n", PetscBools[distributed]);CHKERRQ(ierr);
679*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(v, "DMPlexIsInterpolatedCollective: %s\n", DMPlexInterpolatedFlags[interpolated]);CHKERRQ(ierr);
680*c4762a1bSJed Brown   PetscFunctionReturn(0);
681*c4762a1bSJed Brown }
682*c4762a1bSJed Brown 
683*c4762a1bSJed Brown static PetscErrorCode CreateMeshFromFile(MPI_Comm comm, AppCtx *user, DM *dm, DM *serialDM)
684*c4762a1bSJed Brown {
685*c4762a1bSJed Brown   const char    *filename       = user->filename;
686*c4762a1bSJed Brown   PetscBool      testHeavy      = user->testHeavy;
687*c4762a1bSJed Brown   PetscBool      interpCreate   = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE;
688*c4762a1bSJed Brown   PetscBool      distributed    = PETSC_FALSE;
689*c4762a1bSJed Brown   PetscErrorCode ierr;
690*c4762a1bSJed Brown 
691*c4762a1bSJed Brown   PetscFunctionBegin;
692*c4762a1bSJed Brown   *serialDM = NULL;
693*c4762a1bSJed Brown   if (testHeavy && interpCreate) {ierr = DMPlexSetOrientInterface_Private(NULL, PETSC_FALSE);CHKERRQ(ierr);}
694*c4762a1bSJed Brown   ierr = PetscLogStagePush(stage[0]);CHKERRQ(ierr);
695*c4762a1bSJed Brown   ierr = DMPlexCreateFromFile(comm, filename, interpCreate, dm);CHKERRQ(ierr); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
696*c4762a1bSJed Brown   ierr = PetscLogStagePop();CHKERRQ(ierr);
697*c4762a1bSJed Brown   if (testHeavy && interpCreate) {ierr = DMPlexSetOrientInterface_Private(NULL, PETSC_TRUE);CHKERRQ(ierr);}
698*c4762a1bSJed Brown   ierr = DMPlexIsDistributed(*dm, &distributed);CHKERRQ(ierr);
699*c4762a1bSJed Brown   ierr = PetscPrintf(comm, "DMPlexCreateFromFile produced %s mesh.\n", distributed ? "distributed" : "serial");CHKERRQ(ierr);
700*c4762a1bSJed Brown   if (testHeavy && distributed) {
701*c4762a1bSJed Brown     ierr = PetscOptionsSetValue(NULL, "-dm_plex_hdf5_force_sequential", NULL);CHKERRQ(ierr);
702*c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, filename, interpCreate, serialDM);CHKERRQ(ierr);
703*c4762a1bSJed Brown     ierr = DMPlexIsDistributed(*serialDM, &distributed);CHKERRQ(ierr);
704*c4762a1bSJed Brown     if (distributed) SETERRQ(comm, PETSC_ERR_PLIB, "unable to create a serial DM from file");
705*c4762a1bSJed Brown   }
706*c4762a1bSJed Brown   ierr = DMGetDimension(*dm, &user->dim);CHKERRQ(ierr);
707*c4762a1bSJed Brown   PetscFunctionReturn(0);
708*c4762a1bSJed Brown }
709*c4762a1bSJed Brown 
710*c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
711*c4762a1bSJed Brown {
712*c4762a1bSJed Brown   PetscPartitioner part;
713*c4762a1bSJed Brown   PortableBoundary boundary     = NULL;
714*c4762a1bSJed Brown   DM             serialDM       = NULL;
715*c4762a1bSJed Brown   PetscBool      cellSimplex    = user->cellSimplex;
716*c4762a1bSJed Brown   PetscBool      useGenerator   = user->useGenerator;
717*c4762a1bSJed Brown   PetscBool      interpCreate   = user->interpolate == CREATE ? PETSC_TRUE : PETSC_FALSE;
718*c4762a1bSJed Brown   PetscBool      interpSerial   = user->interpolate == AFTER_CREATE ? PETSC_TRUE : PETSC_FALSE;
719*c4762a1bSJed Brown   PetscBool      interpParallel = user->interpolate == AFTER_DISTRIBUTE ? PETSC_TRUE : PETSC_FALSE;
720*c4762a1bSJed Brown   PetscBool      testHeavy      = user->testHeavy;
721*c4762a1bSJed Brown   PetscMPIInt    rank;
722*c4762a1bSJed Brown   PetscErrorCode ierr;
723*c4762a1bSJed Brown 
724*c4762a1bSJed Brown   PetscFunctionBegin;
725*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
726*c4762a1bSJed Brown   if (user->filename[0]) {
727*c4762a1bSJed Brown     ierr = CreateMeshFromFile(comm, user, dm, &serialDM);CHKERRQ(ierr);
728*c4762a1bSJed Brown   } else if (useGenerator) {
729*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[0]);CHKERRQ(ierr);
730*c4762a1bSJed Brown     ierr = DMPlexCreateBoxMesh(comm, user->dim, cellSimplex, user->faces, NULL, NULL, NULL, interpCreate, dm);CHKERRQ(ierr);
731*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
732*c4762a1bSJed Brown   } else {
733*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[0]);CHKERRQ(ierr);
734*c4762a1bSJed Brown     switch (user->dim) {
735*c4762a1bSJed Brown     case 1:
736*c4762a1bSJed Brown       ierr = CreateMesh_1D(comm, interpCreate, user, dm);CHKERRQ(ierr);
737*c4762a1bSJed Brown       break;
738*c4762a1bSJed Brown     case 2:
739*c4762a1bSJed Brown       if (cellSimplex) {
740*c4762a1bSJed Brown         ierr = CreateSimplex_2D(comm, interpCreate, user, dm);CHKERRQ(ierr);
741*c4762a1bSJed Brown       } else {
742*c4762a1bSJed Brown         ierr = CreateQuad_2D(comm, interpCreate, user, dm);CHKERRQ(ierr);
743*c4762a1bSJed Brown       }
744*c4762a1bSJed Brown       break;
745*c4762a1bSJed Brown     case 3:
746*c4762a1bSJed Brown       if (cellSimplex) {
747*c4762a1bSJed Brown         ierr = CreateSimplex_3D(comm, interpCreate, user, dm);CHKERRQ(ierr);
748*c4762a1bSJed Brown       } else {
749*c4762a1bSJed Brown         ierr = CreateHex_3D(comm, interpCreate, user, dm);CHKERRQ(ierr);
750*c4762a1bSJed Brown       }
751*c4762a1bSJed Brown       break;
752*c4762a1bSJed Brown     default:
753*c4762a1bSJed Brown       SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %D", user->dim);
754*c4762a1bSJed Brown     }
755*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
756*c4762a1bSJed Brown   }
757*c4762a1bSJed Brown   if (user->ncoords % user->dim) SETERRQ2(comm, PETSC_ERR_ARG_OUTOFRANGE, "length of coordinates array %D must be divisable by spatial dimension %D", user->ncoords, user->dim);
758*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Original Mesh");CHKERRQ(ierr);
759*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-orig_dm_view");CHKERRQ(ierr);
760*c4762a1bSJed Brown 
761*c4762a1bSJed Brown   if (interpSerial) {
762*c4762a1bSJed Brown     DM idm;
763*c4762a1bSJed Brown 
764*c4762a1bSJed Brown     if (testHeavy) {ierr = DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE);CHKERRQ(ierr);}
765*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[2]);CHKERRQ(ierr);
766*c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
767*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
768*c4762a1bSJed Brown     if (testHeavy) {ierr = DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE);CHKERRQ(ierr);}
769*c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
770*c4762a1bSJed Brown     *dm = idm;
771*c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) *dm, "Interpolated Mesh");CHKERRQ(ierr);
772*c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-intp_dm_view");CHKERRQ(ierr);
773*c4762a1bSJed Brown   }
774*c4762a1bSJed Brown 
775*c4762a1bSJed Brown   /* Set partitioner options */
776*c4762a1bSJed Brown   ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
777*c4762a1bSJed Brown   if (part) {
778*c4762a1bSJed Brown     ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSIMPLE);CHKERRQ(ierr);
779*c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
780*c4762a1bSJed Brown   }
781*c4762a1bSJed Brown 
782*c4762a1bSJed Brown   if (user->customView) {ierr = CustomView(*dm, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);}
783*c4762a1bSJed Brown   if (testHeavy) {
784*c4762a1bSJed Brown     PetscBool distributed;
785*c4762a1bSJed Brown 
786*c4762a1bSJed Brown     ierr = DMPlexIsDistributed(*dm, &distributed);CHKERRQ(ierr);
787*c4762a1bSJed Brown     if (!serialDM && !distributed) {
788*c4762a1bSJed Brown       serialDM = *dm;
789*c4762a1bSJed Brown       ierr = PetscObjectReference((PetscObject)*dm);CHKERRQ(ierr);
790*c4762a1bSJed Brown     }
791*c4762a1bSJed Brown     if (serialDM) {
792*c4762a1bSJed Brown       ierr = DMPlexGetExpandedBoundary_Private(serialDM, &boundary);CHKERRQ(ierr);
793*c4762a1bSJed Brown     }
794*c4762a1bSJed Brown     if (boundary) {
795*c4762a1bSJed Brown       /* check DM which has been created in parallel and already interpolated */
796*c4762a1bSJed Brown       ierr = DMPlexCheckPointSFHeavy(*dm, boundary);CHKERRQ(ierr);
797*c4762a1bSJed Brown     }
798*c4762a1bSJed Brown     /* Orient interface because it could be deliberately skipped above. It is idempotent. */
799*c4762a1bSJed Brown     ierr = DMPlexOrientInterface_Internal(*dm);CHKERRQ(ierr);
800*c4762a1bSJed Brown   }
801*c4762a1bSJed Brown   if (user->distribute) {
802*c4762a1bSJed Brown     DM               pdm = NULL;
803*c4762a1bSJed Brown 
804*c4762a1bSJed Brown     /* Redistribute mesh over processes using that partitioner */
805*c4762a1bSJed Brown     ierr = PetscLogStagePush(stage[1]);CHKERRQ(ierr);
806*c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
807*c4762a1bSJed Brown     ierr = PetscLogStagePop();CHKERRQ(ierr);
808*c4762a1bSJed Brown     if (pdm) {
809*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
810*c4762a1bSJed Brown       *dm  = pdm;
811*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) *dm, "Redistributed Mesh");CHKERRQ(ierr);
812*c4762a1bSJed Brown       ierr = DMViewFromOptions(*dm, NULL, "-dist_dm_view");CHKERRQ(ierr);
813*c4762a1bSJed Brown     }
814*c4762a1bSJed Brown 
815*c4762a1bSJed Brown     if (interpParallel) {
816*c4762a1bSJed Brown       DM idm;
817*c4762a1bSJed Brown 
818*c4762a1bSJed Brown       if (testHeavy) {ierr = DMPlexSetOrientInterface_Private(*dm, PETSC_FALSE);CHKERRQ(ierr);}
819*c4762a1bSJed Brown       ierr = PetscLogStagePush(stage[2]);CHKERRQ(ierr);
820*c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexCheckPointSFHeavy() */
821*c4762a1bSJed Brown       ierr = PetscLogStagePop();CHKERRQ(ierr);
822*c4762a1bSJed Brown       if (testHeavy) {ierr = DMPlexSetOrientInterface_Private(*dm, PETSC_TRUE);CHKERRQ(ierr);}
823*c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
824*c4762a1bSJed Brown       *dm = idm;
825*c4762a1bSJed Brown       ierr = PetscObjectSetName((PetscObject) *dm, "Interpolated Redistributed Mesh");CHKERRQ(ierr);
826*c4762a1bSJed Brown       ierr = DMViewFromOptions(*dm, NULL, "-intp_dm_view");CHKERRQ(ierr);
827*c4762a1bSJed Brown     }
828*c4762a1bSJed Brown   }
829*c4762a1bSJed Brown   if (testHeavy) {
830*c4762a1bSJed Brown     if (boundary) {
831*c4762a1bSJed Brown       ierr = DMPlexCheckPointSFHeavy(*dm, boundary);CHKERRQ(ierr);
832*c4762a1bSJed Brown     }
833*c4762a1bSJed Brown     /* Orient interface because it could be deliberately skipped above. It is idempotent. */
834*c4762a1bSJed Brown     ierr = DMPlexOrientInterface_Internal(*dm);CHKERRQ(ierr);
835*c4762a1bSJed Brown   }
836*c4762a1bSJed Brown 
837*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Parallel Mesh");CHKERRQ(ierr);
838*c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
839*c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
840*c4762a1bSJed Brown 
841*c4762a1bSJed Brown   if (user->customView) {ierr = CustomView(*dm, PETSC_VIEWER_STDOUT_(comm));CHKERRQ(ierr);}
842*c4762a1bSJed Brown   ierr = DMDestroy(&serialDM);CHKERRQ(ierr);
843*c4762a1bSJed Brown   ierr = PortableBoundaryDestroy(&boundary);CHKERRQ(ierr);
844*c4762a1bSJed Brown   PetscFunctionReturn(0);
845*c4762a1bSJed Brown }
846*c4762a1bSJed Brown 
847*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscErrorCode coord2str(char buf[], size_t len, PetscInt dim, PetscReal *coords, PetscReal tol)
848*c4762a1bSJed Brown {
849*c4762a1bSJed Brown   PetscErrorCode ierr;
850*c4762a1bSJed Brown 
851*c4762a1bSJed Brown   PetscFunctionBegin;
852*c4762a1bSJed Brown   if (dim > 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "dim must be less than or equal 3");
853*c4762a1bSJed Brown   if (tol >= 1e-3) {
854*c4762a1bSJed Brown     switch (dim) {
855*c4762a1bSJed Brown       case 1: ierr = PetscSNPrintf(buf,len,"(%12.3f)",(double)coords[0]);CHKERRQ(ierr);
856*c4762a1bSJed Brown       case 2: ierr = PetscSNPrintf(buf,len,"(%12.3f, %12.3f)",(double)coords[0],(double)coords[1]);CHKERRQ(ierr);
857*c4762a1bSJed Brown       case 3: ierr = PetscSNPrintf(buf,len,"(%12.3f, %12.3f, %12.3f)",(double)coords[0],(double)coords[1],(double)coords[2]);CHKERRQ(ierr);
858*c4762a1bSJed Brown     }
859*c4762a1bSJed Brown   } else {
860*c4762a1bSJed Brown     switch (dim) {
861*c4762a1bSJed Brown       case 1: ierr = PetscSNPrintf(buf,len,"(%12.6f)",(double)coords[0]);CHKERRQ(ierr);
862*c4762a1bSJed Brown       case 2: ierr = PetscSNPrintf(buf,len,"(%12.6f, %12.6f)",(double)coords[0],(double)coords[1]);CHKERRQ(ierr);
863*c4762a1bSJed Brown       case 3: ierr = PetscSNPrintf(buf,len,"(%12.6f, %12.6f, %12.6f)",(double)coords[0],(double)coords[1],(double)coords[2]);CHKERRQ(ierr);
864*c4762a1bSJed Brown     }
865*c4762a1bSJed Brown   }
866*c4762a1bSJed Brown   PetscFunctionReturn(0);
867*c4762a1bSJed Brown }
868*c4762a1bSJed Brown 
869*c4762a1bSJed Brown static PetscErrorCode ViewVerticesFromCoords(DM dm, PetscInt npoints, PetscReal coords[], PetscReal tol, PetscViewer viewer)
870*c4762a1bSJed Brown {
871*c4762a1bSJed Brown   PetscInt       dim, i;
872*c4762a1bSJed Brown   PetscInt       *points;
873*c4762a1bSJed Brown   char           coordstr[128];
874*c4762a1bSJed Brown   MPI_Comm       comm;
875*c4762a1bSJed Brown   PetscMPIInt    rank;
876*c4762a1bSJed Brown   PetscErrorCode ierr;
877*c4762a1bSJed Brown 
878*c4762a1bSJed Brown   PetscFunctionBegin;
879*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
880*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
881*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
882*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr);
883*c4762a1bSJed Brown   ierr = PetscMalloc1(npoints, &points);CHKERRQ(ierr);
884*c4762a1bSJed Brown   ierr = DMPlexFindVertices(dm, npoints, coords, tol, points);CHKERRQ(ierr);
885*c4762a1bSJed Brown   for (i=0; i < npoints; i++) {
886*c4762a1bSJed Brown     ierr = coord2str(coordstr, sizeof(coordstr), dim, &coords[i*dim], tol);CHKERRQ(ierr);
887*c4762a1bSJed Brown     if (!rank && i) {ierr = PetscViewerASCIISynchronizedPrintf(viewer, "-----\n");CHKERRQ(ierr);}
888*c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %s --> points[%D] = %D\n", rank, coordstr, i, points[i]);CHKERRQ(ierr);
889*c4762a1bSJed Brown     ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
890*c4762a1bSJed Brown   }
891*c4762a1bSJed Brown   ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr);
892*c4762a1bSJed Brown   ierr = PetscFree(points);CHKERRQ(ierr);
893*c4762a1bSJed Brown   PetscFunctionReturn(0);
894*c4762a1bSJed Brown }
895*c4762a1bSJed Brown 
896*c4762a1bSJed Brown static PetscErrorCode TestExpandPoints(DM dm, AppCtx *user)
897*c4762a1bSJed Brown {
898*c4762a1bSJed Brown   IS                is;
899*c4762a1bSJed Brown   PetscSection      *sects;
900*c4762a1bSJed Brown   IS                *iss;
901*c4762a1bSJed Brown   PetscInt          d,depth;
902*c4762a1bSJed Brown   PetscMPIInt       rank;
903*c4762a1bSJed Brown   PetscErrorCode    ierr;
904*c4762a1bSJed Brown   PetscViewer       viewer=PETSC_VIEWER_STDOUT_WORLD, sviewer;
905*c4762a1bSJed Brown 
906*c4762a1bSJed Brown   PetscFunctionBegin;
907*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
908*c4762a1bSJed Brown   if (user->testExpandPointsEmpty && !rank) {
909*c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &is);CHKERRQ(ierr);
910*c4762a1bSJed Brown   } else {
911*c4762a1bSJed Brown     ierr = ISCreateGeneral(PETSC_COMM_SELF, user->nPointsToExpand, user->pointsToExpand, PETSC_USE_POINTER, &is);CHKERRQ(ierr);
912*c4762a1bSJed Brown   }
913*c4762a1bSJed Brown   ierr = DMPlexGetConeRecursive(dm, is, &depth, &iss, &sects);CHKERRQ(ierr);
914*c4762a1bSJed Brown   ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
915*c4762a1bSJed Brown   ierr = PetscViewerASCIIPrintf(sviewer, "[%d] ==========================\n",rank);CHKERRQ(ierr);
916*c4762a1bSJed Brown   for (d=depth-1; d>=0; d--) {
917*c4762a1bSJed Brown     IS          checkIS;
918*c4762a1bSJed Brown     PetscBool   flg;
919*c4762a1bSJed Brown 
920*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(sviewer, "depth %D ---------------\n",d);CHKERRQ(ierr);
921*c4762a1bSJed Brown     ierr = PetscSectionView(sects[d], sviewer);CHKERRQ(ierr);
922*c4762a1bSJed Brown     ierr = ISView(iss[d], sviewer);CHKERRQ(ierr);
923*c4762a1bSJed Brown     /* check reverse operation */
924*c4762a1bSJed Brown     if (d < depth-1) {
925*c4762a1bSJed Brown       ierr = DMPlexExpandedConesToFaces_Private(dm, iss[d], sects[d], &checkIS);CHKERRQ(ierr);
926*c4762a1bSJed Brown       ierr = ISEqualUnsorted(checkIS, iss[d+1], &flg);CHKERRQ(ierr);
927*c4762a1bSJed Brown       if (!flg) SETERRQ(PetscObjectComm((PetscObject) checkIS), PETSC_ERR_PLIB, "DMPlexExpandedConesToFaces_Private produced wrong IS");
928*c4762a1bSJed Brown       ierr = ISDestroy(&checkIS);CHKERRQ(ierr);
929*c4762a1bSJed Brown     }
930*c4762a1bSJed Brown   }
931*c4762a1bSJed Brown   ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
932*c4762a1bSJed Brown   ierr = PetscViewerFlush(viewer);CHKERRQ(ierr);
933*c4762a1bSJed Brown   ierr = DMPlexRestoreConeRecursive(dm, is, &depth, &iss, &sects);CHKERRQ(ierr);
934*c4762a1bSJed Brown   ierr = ISDestroy(&is);CHKERRQ(ierr);
935*c4762a1bSJed Brown   PetscFunctionReturn(0);
936*c4762a1bSJed Brown }
937*c4762a1bSJed Brown 
938*c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedConesToFaces_Private(DM dm, IS is, PetscSection section, IS *newis)
939*c4762a1bSJed Brown {
940*c4762a1bSJed Brown   PetscInt          n,n1,ncone,numCoveredPoints,o,p,q,start,end;
941*c4762a1bSJed Brown   const PetscInt    *coveredPoints;
942*c4762a1bSJed Brown   const PetscInt    *arr, *cone;
943*c4762a1bSJed Brown   PetscInt          *newarr;
944*c4762a1bSJed Brown   PetscErrorCode    ierr;
945*c4762a1bSJed Brown 
946*c4762a1bSJed Brown   PetscFunctionBegin;
947*c4762a1bSJed Brown   ierr = ISGetLocalSize(is, &n);CHKERRQ(ierr);
948*c4762a1bSJed Brown   ierr = PetscSectionGetStorageSize(section, &n1);CHKERRQ(ierr);
949*c4762a1bSJed Brown   ierr = PetscSectionGetChart(section, &start, &end);CHKERRQ(ierr);
950*c4762a1bSJed Brown   if (n != n1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "IS size = %D != %D = section storage size\n", n, n1);
951*c4762a1bSJed Brown   ierr = ISGetIndices(is, &arr);CHKERRQ(ierr);
952*c4762a1bSJed Brown   ierr = PetscMalloc1(end-start, &newarr);CHKERRQ(ierr);
953*c4762a1bSJed Brown   for (q=start; q<end; q++) {
954*c4762a1bSJed Brown     ierr = PetscSectionGetDof(section, q, &ncone);CHKERRQ(ierr);
955*c4762a1bSJed Brown     ierr = PetscSectionGetOffset(section, q, &o);CHKERRQ(ierr);
956*c4762a1bSJed Brown     cone = &arr[o];
957*c4762a1bSJed Brown     if (ncone == 1) {
958*c4762a1bSJed Brown       numCoveredPoints = 1;
959*c4762a1bSJed Brown       p = cone[0];
960*c4762a1bSJed Brown     } else {
961*c4762a1bSJed Brown       PetscInt i;
962*c4762a1bSJed Brown       p = PETSC_MAX_INT;
963*c4762a1bSJed Brown       for (i=0; i<ncone; i++) if (cone[i] < 0) {p = -1; break;}
964*c4762a1bSJed Brown       if (p >= 0) {
965*c4762a1bSJed Brown         ierr = DMPlexGetJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
966*c4762a1bSJed Brown         if (numCoveredPoints > 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "more than one covered points for section point %D",q);
967*c4762a1bSJed Brown         if (numCoveredPoints) p = coveredPoints[0];
968*c4762a1bSJed Brown         else                  p = -2;
969*c4762a1bSJed Brown         ierr = DMPlexRestoreJoin(dm, ncone, cone, &numCoveredPoints, &coveredPoints);CHKERRQ(ierr);
970*c4762a1bSJed Brown       }
971*c4762a1bSJed Brown     }
972*c4762a1bSJed Brown     newarr[q-start] = p;
973*c4762a1bSJed Brown   }
974*c4762a1bSJed Brown   ierr = ISRestoreIndices(is, &arr);CHKERRQ(ierr);
975*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF, end-start, newarr, PETSC_OWN_POINTER, newis);CHKERRQ(ierr);
976*c4762a1bSJed Brown   PetscFunctionReturn(0);
977*c4762a1bSJed Brown }
978*c4762a1bSJed Brown 
979*c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedVerticesToFaces_Private(DM dm, IS boundary_expanded_is, PetscInt depth, PetscSection sections[], IS *boundary_is)
980*c4762a1bSJed Brown {
981*c4762a1bSJed Brown   PetscInt          d;
982*c4762a1bSJed Brown   IS                is,newis;
983*c4762a1bSJed Brown   PetscErrorCode    ierr;
984*c4762a1bSJed Brown 
985*c4762a1bSJed Brown   PetscFunctionBegin;
986*c4762a1bSJed Brown   is = boundary_expanded_is;
987*c4762a1bSJed Brown   ierr = PetscObjectReference((PetscObject)is);CHKERRQ(ierr);
988*c4762a1bSJed Brown   for (d = 0; d < depth-1; ++d) {
989*c4762a1bSJed Brown     ierr = DMPlexExpandedConesToFaces_Private(dm, is, sections[d], &newis);CHKERRQ(ierr);
990*c4762a1bSJed Brown     ierr = ISDestroy(&is);CHKERRQ(ierr);
991*c4762a1bSJed Brown     is = newis;
992*c4762a1bSJed Brown   }
993*c4762a1bSJed Brown   *boundary_is = is;
994*c4762a1bSJed Brown   PetscFunctionReturn(0);
995*c4762a1bSJed Brown }
996*c4762a1bSJed Brown 
997*c4762a1bSJed Brown #define CHKERRQI(incall,ierr) if (ierr) {incall = PETSC_FALSE; }
998*c4762a1bSJed Brown 
999*c4762a1bSJed Brown static PetscErrorCode DMLabelViewFromOptionsOnComm_Private(DMLabel label, const char optionname[], MPI_Comm comm)
1000*c4762a1bSJed Brown {
1001*c4762a1bSJed Brown   PetscErrorCode    ierr;
1002*c4762a1bSJed Brown   PetscViewer       viewer;
1003*c4762a1bSJed Brown   PetscBool         flg;
1004*c4762a1bSJed Brown   static PetscBool  incall = PETSC_FALSE;
1005*c4762a1bSJed Brown   PetscViewerFormat format;
1006*c4762a1bSJed Brown 
1007*c4762a1bSJed Brown   PetscFunctionBegin;
1008*c4762a1bSJed Brown   if (incall) PetscFunctionReturn(0);
1009*c4762a1bSJed Brown   incall = PETSC_TRUE;
1010*c4762a1bSJed Brown   ierr   = PetscOptionsGetViewer(comm,((PetscObject)label)->options,((PetscObject)label)->prefix,optionname,&viewer,&format,&flg);CHKERRQI(incall,ierr);
1011*c4762a1bSJed Brown   if (flg) {
1012*c4762a1bSJed Brown     ierr = PetscViewerPushFormat(viewer,format);CHKERRQI(incall,ierr);
1013*c4762a1bSJed Brown     ierr = DMLabelView(label, viewer);CHKERRQI(incall,ierr);
1014*c4762a1bSJed Brown     ierr = PetscViewerPopFormat(viewer);CHKERRQI(incall,ierr);
1015*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&viewer);CHKERRQI(incall,ierr);
1016*c4762a1bSJed Brown   }
1017*c4762a1bSJed Brown   incall = PETSC_FALSE;
1018*c4762a1bSJed Brown   PetscFunctionReturn(0);
1019*c4762a1bSJed Brown }
1020*c4762a1bSJed Brown 
1021*c4762a1bSJed Brown /* TODO: this is hotfixing DMLabelGetStratumIS() - it should be fixed systematically instead */
1022*c4762a1bSJed Brown PETSC_STATIC_INLINE PetscErrorCode DMLabelGetStratumISOnComm_Private(DMLabel label, PetscInt value, MPI_Comm comm, IS *is)
1023*c4762a1bSJed Brown {
1024*c4762a1bSJed Brown   IS                tmpis;
1025*c4762a1bSJed Brown   PetscErrorCode    ierr;
1026*c4762a1bSJed Brown 
1027*c4762a1bSJed Brown   PetscFunctionBegin;
1028*c4762a1bSJed Brown   ierr = DMLabelGetStratumIS(label, value, &tmpis);CHKERRQ(ierr);
1029*c4762a1bSJed Brown   if (!tmpis) {ierr = ISCreateGeneral(PETSC_COMM_SELF, 0, NULL, PETSC_USE_POINTER, &tmpis);CHKERRQ(ierr);}
1030*c4762a1bSJed Brown   ierr = ISOnComm(tmpis, comm, PETSC_COPY_VALUES, is);CHKERRQ(ierr);
1031*c4762a1bSJed Brown   ierr = ISDestroy(&tmpis);CHKERRQ(ierr);
1032*c4762a1bSJed Brown   PetscFunctionReturn(0);
1033*c4762a1bSJed Brown }
1034*c4762a1bSJed Brown 
1035*c4762a1bSJed Brown /* currently only for simple PetscSection without fields or constraints */
1036*c4762a1bSJed Brown static PetscErrorCode PetscSectionReplicate_Private(MPI_Comm comm, PetscMPIInt rootrank, PetscSection sec0, PetscSection *secout)
1037*c4762a1bSJed Brown {
1038*c4762a1bSJed Brown   PetscSection      sec;
1039*c4762a1bSJed Brown   PetscInt          chart[2], p;
1040*c4762a1bSJed Brown   PetscInt          *dofarr;
1041*c4762a1bSJed Brown   PetscMPIInt       rank;
1042*c4762a1bSJed Brown   PetscErrorCode    ierr;
1043*c4762a1bSJed Brown 
1044*c4762a1bSJed Brown   PetscFunctionBegin;
1045*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1046*c4762a1bSJed Brown   if (rank == rootrank) {
1047*c4762a1bSJed Brown     ierr = PetscSectionGetChart(sec0, &chart[0], &chart[1]);CHKERRQ(ierr);
1048*c4762a1bSJed Brown   }
1049*c4762a1bSJed Brown   ierr = MPI_Bcast(chart, 2, MPIU_INT, rootrank, comm);CHKERRQ(ierr);
1050*c4762a1bSJed Brown   ierr = PetscMalloc1(chart[1]-chart[0], &dofarr);CHKERRQ(ierr);
1051*c4762a1bSJed Brown   if (rank == rootrank) {
1052*c4762a1bSJed Brown     for (p = chart[0]; p < chart[1]; p++) {
1053*c4762a1bSJed Brown       ierr = PetscSectionGetDof(sec0, p, &dofarr[p-chart[0]]);CHKERRQ(ierr);
1054*c4762a1bSJed Brown     }
1055*c4762a1bSJed Brown   }
1056*c4762a1bSJed Brown   ierr = MPI_Bcast(dofarr, chart[1]-chart[0], MPIU_INT, rootrank, comm);CHKERRQ(ierr);
1057*c4762a1bSJed Brown   ierr = PetscSectionCreate(comm, &sec);CHKERRQ(ierr);
1058*c4762a1bSJed Brown   ierr = PetscSectionSetChart(sec, chart[0], chart[1]);CHKERRQ(ierr);
1059*c4762a1bSJed Brown   for (p = chart[0]; p < chart[1]; p++) {
1060*c4762a1bSJed Brown     ierr = PetscSectionSetDof(sec, p, dofarr[p-chart[0]]);CHKERRQ(ierr);
1061*c4762a1bSJed Brown   }
1062*c4762a1bSJed Brown   ierr = PetscSectionSetUp(sec);CHKERRQ(ierr);
1063*c4762a1bSJed Brown   ierr = PetscFree(dofarr);CHKERRQ(ierr);
1064*c4762a1bSJed Brown   *secout = sec;
1065*c4762a1bSJed Brown   PetscFunctionReturn(0);
1066*c4762a1bSJed Brown }
1067*c4762a1bSJed Brown 
1068*c4762a1bSJed Brown static PetscErrorCode VecToPetscReal_Private(Vec vec, PetscReal *rvals[])
1069*c4762a1bSJed Brown {
1070*c4762a1bSJed Brown   PetscInt          n;
1071*c4762a1bSJed Brown   const PetscScalar *svals;
1072*c4762a1bSJed Brown   PetscErrorCode    ierr;
1073*c4762a1bSJed Brown 
1074*c4762a1bSJed Brown   PetscFunctionBegin;
1075*c4762a1bSJed Brown   ierr = VecGetLocalSize(vec, &n);CHKERRQ(ierr);
1076*c4762a1bSJed Brown   ierr = VecGetArrayRead(vec, &svals);CHKERRQ(ierr);
1077*c4762a1bSJed Brown   ierr = PetscMalloc1(n, rvals);CHKERRQ(ierr);
1078*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
1079*c4762a1bSJed Brown   {
1080*c4762a1bSJed Brown     PetscInt i;
1081*c4762a1bSJed Brown     for (i=0; i<n; i++) (*rvals)[i] = PetscRealPart(svals[i]);
1082*c4762a1bSJed Brown   }
1083*c4762a1bSJed Brown #else
1084*c4762a1bSJed Brown   ierr = PetscMemcpy(*rvals, svals, n*sizeof(PetscReal));
1085*c4762a1bSJed Brown #endif
1086*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(vec, &svals);CHKERRQ(ierr);
1087*c4762a1bSJed Brown   PetscFunctionReturn(0);
1088*c4762a1bSJed Brown }
1089*c4762a1bSJed Brown 
1090*c4762a1bSJed Brown static PetscErrorCode DMPlexExpandedVerticesCoordinatesToFaces_Private(DM ipdm, PortableBoundary bnd, IS *face_is)
1091*c4762a1bSJed Brown {
1092*c4762a1bSJed Brown   PetscInt            dim, ncoords, npoints;
1093*c4762a1bSJed Brown   PetscReal           *rcoords;
1094*c4762a1bSJed Brown   PetscInt            *points;
1095*c4762a1bSJed Brown   IS                  faces_expanded_is;
1096*c4762a1bSJed Brown   PetscErrorCode      ierr;
1097*c4762a1bSJed Brown 
1098*c4762a1bSJed Brown   PetscFunctionBegin;
1099*c4762a1bSJed Brown   ierr = DMGetDimension(ipdm, &dim);CHKERRQ(ierr);
1100*c4762a1bSJed Brown   ierr = VecGetLocalSize(bnd->coordinates, &ncoords);CHKERRQ(ierr);
1101*c4762a1bSJed Brown   ierr = VecToPetscReal_Private(bnd->coordinates, &rcoords);CHKERRQ(ierr);
1102*c4762a1bSJed Brown   npoints = ncoords / dim;
1103*c4762a1bSJed Brown   ierr = PetscMalloc1(npoints, &points);CHKERRQ(ierr);
1104*c4762a1bSJed Brown   ierr = DMPlexFindVertices(ipdm, npoints, rcoords, 0.0, points);CHKERRQ(ierr);
1105*c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF, npoints, points, PETSC_OWN_POINTER, &faces_expanded_is);CHKERRQ(ierr);
1106*c4762a1bSJed Brown   ierr = DMPlexExpandedVerticesToFaces_Private(ipdm, faces_expanded_is, bnd->depth, bnd->sections, face_is);CHKERRQ(ierr);
1107*c4762a1bSJed Brown   ierr = PetscFree(rcoords);CHKERRQ(ierr);
1108*c4762a1bSJed Brown   ierr = ISDestroy(&faces_expanded_is);CHKERRQ(ierr);
1109*c4762a1bSJed Brown   PetscFunctionReturn(0);
1110*c4762a1bSJed Brown }
1111*c4762a1bSJed Brown 
1112*c4762a1bSJed Brown /* hack disabling DMPlexOrientInterface() call in DMPlexInterpolate() via -dm_plex_interpolate_orient_interfaces option */
1113*c4762a1bSJed Brown static PetscErrorCode DMPlexSetOrientInterface_Private(DM dm, PetscBool enable)
1114*c4762a1bSJed Brown {
1115*c4762a1bSJed Brown   PetscOptions      options = NULL;
1116*c4762a1bSJed Brown   const char        *prefix = NULL;
1117*c4762a1bSJed Brown   const char        opt[] = "-dm_plex_interpolate_orient_interfaces";
1118*c4762a1bSJed Brown   char              prefix_opt[512];
1119*c4762a1bSJed Brown   PetscBool         flg, set;
1120*c4762a1bSJed Brown   static PetscBool  wasSetTrue = PETSC_FALSE;
1121*c4762a1bSJed Brown   PetscErrorCode    ierr;
1122*c4762a1bSJed Brown 
1123*c4762a1bSJed Brown   PetscFunctionBegin;
1124*c4762a1bSJed Brown   if (dm) {
1125*c4762a1bSJed Brown     ierr = PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix);CHKERRQ(ierr);
1126*c4762a1bSJed Brown     options = ((PetscObject)dm)->options;
1127*c4762a1bSJed Brown   }
1128*c4762a1bSJed Brown   ierr = PetscStrcpy(prefix_opt, "-");CHKERRQ(ierr);
1129*c4762a1bSJed Brown   ierr = PetscStrlcat(prefix_opt, prefix, sizeof(prefix_opt));CHKERRQ(ierr);
1130*c4762a1bSJed Brown   ierr = PetscStrlcat(prefix_opt, &opt[1], sizeof(prefix_opt));CHKERRQ(ierr);
1131*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(options, prefix, opt, &flg, &set);CHKERRQ(ierr);
1132*c4762a1bSJed Brown   if (!enable) {
1133*c4762a1bSJed Brown     if (set && flg) wasSetTrue = PETSC_TRUE;
1134*c4762a1bSJed Brown     ierr = PetscOptionsSetValue(options, prefix_opt, "0");CHKERRQ(ierr);
1135*c4762a1bSJed Brown   } else if (set && !flg) {
1136*c4762a1bSJed Brown     if (wasSetTrue) {
1137*c4762a1bSJed Brown       ierr = PetscOptionsSetValue(options, prefix_opt, "1");CHKERRQ(ierr);
1138*c4762a1bSJed Brown     } else {
1139*c4762a1bSJed Brown       /* default is PETSC_TRUE */
1140*c4762a1bSJed Brown       ierr = PetscOptionsClearValue(options, prefix_opt);CHKERRQ(ierr);
1141*c4762a1bSJed Brown     }
1142*c4762a1bSJed Brown     wasSetTrue = PETSC_FALSE;
1143*c4762a1bSJed Brown   }
1144*c4762a1bSJed Brown #if defined(PETSC_USE_DEBUG)
1145*c4762a1bSJed Brown   {
1146*c4762a1bSJed Brown     ierr = PetscOptionsGetBool(options, prefix, opt, &flg, &set);CHKERRQ(ierr);
1147*c4762a1bSJed Brown     if (PetscUnlikely(set && flg != enable)) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "PetscOptionsSetValue did not have the desired effect");
1148*c4762a1bSJed Brown   }
1149*c4762a1bSJed Brown #endif
1150*c4762a1bSJed Brown   PetscFunctionReturn(0);
1151*c4762a1bSJed Brown }
1152*c4762a1bSJed Brown 
1153*c4762a1bSJed Brown /* get coordinate description of the whole-domain boundary */
1154*c4762a1bSJed Brown static PetscErrorCode DMPlexGetExpandedBoundary_Private(DM dm, PortableBoundary *boundary)
1155*c4762a1bSJed Brown {
1156*c4762a1bSJed Brown   PortableBoundary       bnd0, bnd;
1157*c4762a1bSJed Brown   MPI_Comm               comm;
1158*c4762a1bSJed Brown   DM                     idm;
1159*c4762a1bSJed Brown   DMLabel                label;
1160*c4762a1bSJed Brown   PetscInt               d;
1161*c4762a1bSJed Brown   const char             boundaryName[] = "DMPlexDistributeInterpolateMarkInterface_boundary";
1162*c4762a1bSJed Brown   IS                     boundary_is;
1163*c4762a1bSJed Brown   IS                     *boundary_expanded_iss;
1164*c4762a1bSJed Brown   PetscMPIInt            rootrank = 0;
1165*c4762a1bSJed Brown   PetscMPIInt            rank, size;
1166*c4762a1bSJed Brown   PetscInt               value = 1;
1167*c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1168*c4762a1bSJed Brown   PetscBool              flg;
1169*c4762a1bSJed Brown   PetscErrorCode         ierr;
1170*c4762a1bSJed Brown 
1171*c4762a1bSJed Brown   PetscFunctionBegin;
1172*c4762a1bSJed Brown   ierr = PetscNew(&bnd);CHKERRQ(ierr);
1173*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1174*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1175*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1176*c4762a1bSJed Brown   ierr = DMPlexIsDistributed(dm, &flg);CHKERRQ(ierr);
1177*c4762a1bSJed Brown   if (flg) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "serial DM (all points on one rank) needed");
1178*c4762a1bSJed Brown 
1179*c4762a1bSJed Brown   /* interpolate serial DM if not yet interpolated */
1180*c4762a1bSJed Brown   ierr = DMPlexIsInterpolatedCollective(dm, &intp);CHKERRQ(ierr);
1181*c4762a1bSJed Brown   if (intp == DMPLEX_INTERPOLATED_FULL) {
1182*c4762a1bSJed Brown     idm = dm;
1183*c4762a1bSJed Brown     ierr = PetscObjectReference((PetscObject)dm);CHKERRQ(ierr);
1184*c4762a1bSJed Brown   } else {
1185*c4762a1bSJed Brown     ierr = DMPlexInterpolate(dm, &idm);CHKERRQ(ierr);
1186*c4762a1bSJed Brown     ierr = DMViewFromOptions(idm, NULL, "-idm_view");CHKERRQ(ierr);
1187*c4762a1bSJed Brown   }
1188*c4762a1bSJed Brown 
1189*c4762a1bSJed Brown   /* mark whole-domain boundary of the serial DM */
1190*c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, boundaryName, &label);CHKERRQ(ierr);
1191*c4762a1bSJed Brown   ierr = DMAddLabel(idm, label);CHKERRQ(ierr);
1192*c4762a1bSJed Brown   ierr = DMPlexMarkBoundaryFaces(idm, value, label);CHKERRQ(ierr);
1193*c4762a1bSJed Brown   ierr = DMLabelViewFromOptionsOnComm_Private(label, "-idm_boundary_view", comm);CHKERRQ(ierr);
1194*c4762a1bSJed Brown   ierr = DMLabelGetStratumIS(label, value, &boundary_is);CHKERRQ(ierr);
1195*c4762a1bSJed Brown 
1196*c4762a1bSJed Brown   /* translate to coordinates */
1197*c4762a1bSJed Brown   ierr = PetscNew(&bnd0);CHKERRQ(ierr);
1198*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocalSetUp(idm);CHKERRQ(ierr);
1199*c4762a1bSJed Brown   if (rank == rootrank) {
1200*c4762a1bSJed Brown     ierr = DMPlexGetConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections);CHKERRQ(ierr);
1201*c4762a1bSJed Brown     ierr = DMGetCoordinatesLocalTuple(dm, boundary_expanded_iss[0], NULL, &bnd0->coordinates);CHKERRQ(ierr);
1202*c4762a1bSJed Brown     /* self-check */
1203*c4762a1bSJed Brown     {
1204*c4762a1bSJed Brown       IS is0;
1205*c4762a1bSJed Brown       ierr = DMPlexExpandedVerticesCoordinatesToFaces_Private(idm, bnd0, &is0);CHKERRQ(ierr);
1206*c4762a1bSJed Brown       ierr = ISEqual(is0, boundary_is, &flg);CHKERRQ(ierr);
1207*c4762a1bSJed Brown       if (!flg) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "DMPlexExpandedVerticesCoordinatesToFaces_Private produced a wrong IS");
1208*c4762a1bSJed Brown       ierr = ISDestroy(&is0);CHKERRQ(ierr);
1209*c4762a1bSJed Brown     }
1210*c4762a1bSJed Brown   } else {
1211*c4762a1bSJed Brown     ierr = VecCreateSeq(PETSC_COMM_SELF, 0, &bnd0->coordinates);CHKERRQ(ierr);
1212*c4762a1bSJed Brown   }
1213*c4762a1bSJed Brown 
1214*c4762a1bSJed Brown   {
1215*c4762a1bSJed Brown     Vec         tmp;
1216*c4762a1bSJed Brown     VecScatter  sc;
1217*c4762a1bSJed Brown     IS          xis;
1218*c4762a1bSJed Brown     PetscInt    n;
1219*c4762a1bSJed Brown 
1220*c4762a1bSJed Brown     /* just convert seq vectors to mpi vector */
1221*c4762a1bSJed Brown     ierr = VecGetLocalSize(bnd0->coordinates, &n);CHKERRQ(ierr);
1222*c4762a1bSJed Brown     ierr = MPI_Bcast(&n, 1, MPIU_INT, rootrank, comm);CHKERRQ(ierr);
1223*c4762a1bSJed Brown     if (rank == rootrank) {
1224*c4762a1bSJed Brown       ierr = VecCreateMPI(comm, n, n, &tmp);CHKERRQ(ierr);
1225*c4762a1bSJed Brown     } else {
1226*c4762a1bSJed Brown       ierr = VecCreateMPI(comm, 0, n, &tmp);CHKERRQ(ierr);
1227*c4762a1bSJed Brown     }
1228*c4762a1bSJed Brown     ierr = VecCopy(bnd0->coordinates, tmp);CHKERRQ(ierr);
1229*c4762a1bSJed Brown     ierr = VecDestroy(&bnd0->coordinates);CHKERRQ(ierr);
1230*c4762a1bSJed Brown     bnd0->coordinates = tmp;
1231*c4762a1bSJed Brown 
1232*c4762a1bSJed Brown     /* replicate coordinates from root rank to all ranks */
1233*c4762a1bSJed Brown     ierr = VecCreateMPI(comm, n, n*size, &bnd->coordinates);CHKERRQ(ierr);
1234*c4762a1bSJed Brown     ierr = ISCreateStride(comm, n, 0, 1, &xis);CHKERRQ(ierr);
1235*c4762a1bSJed Brown     ierr = VecScatterCreate(bnd0->coordinates, xis, bnd->coordinates, NULL, &sc);CHKERRQ(ierr);
1236*c4762a1bSJed Brown     ierr = VecScatterBegin(sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
1237*c4762a1bSJed Brown     ierr = VecScatterEnd(  sc, bnd0->coordinates, bnd->coordinates, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr);
1238*c4762a1bSJed Brown     ierr = VecScatterDestroy(&sc);CHKERRQ(ierr);
1239*c4762a1bSJed Brown     ierr = ISDestroy(&xis);CHKERRQ(ierr);
1240*c4762a1bSJed Brown   }
1241*c4762a1bSJed Brown   bnd->depth = bnd0->depth;
1242*c4762a1bSJed Brown   ierr = MPI_Bcast(&bnd->depth, 1, MPIU_INT, rootrank, comm);CHKERRQ(ierr);
1243*c4762a1bSJed Brown   ierr = PetscMalloc1(bnd->depth, &bnd->sections);CHKERRQ(ierr);
1244*c4762a1bSJed Brown   for (d=0; d<bnd->depth; d++) {
1245*c4762a1bSJed Brown     ierr = PetscSectionReplicate_Private(comm, rootrank, (rank == rootrank) ? bnd0->sections[d] : NULL, &bnd->sections[d]);CHKERRQ(ierr);
1246*c4762a1bSJed Brown   }
1247*c4762a1bSJed Brown 
1248*c4762a1bSJed Brown   if (rank == rootrank) {
1249*c4762a1bSJed Brown     ierr = DMPlexRestoreConeRecursive(idm, boundary_is, &bnd0->depth, &boundary_expanded_iss, &bnd0->sections);CHKERRQ(ierr);
1250*c4762a1bSJed Brown   }
1251*c4762a1bSJed Brown   ierr = PortableBoundaryDestroy(&bnd0);CHKERRQ(ierr);
1252*c4762a1bSJed Brown   ierr = DMRemoveLabelBySelf(idm, &label, PETSC_TRUE);CHKERRQ(ierr);
1253*c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
1254*c4762a1bSJed Brown   ierr = ISDestroy(&boundary_is);CHKERRQ(ierr);
1255*c4762a1bSJed Brown   ierr = DMDestroy(&idm);CHKERRQ(ierr);
1256*c4762a1bSJed Brown   *boundary = bnd;
1257*c4762a1bSJed Brown   PetscFunctionReturn(0);
1258*c4762a1bSJed Brown }
1259*c4762a1bSJed Brown 
1260*c4762a1bSJed Brown /* get faces of inter-partition interface */
1261*c4762a1bSJed Brown static PetscErrorCode DMPlexGetInterfaceFaces_Private(DM ipdm, IS boundary_faces_is, IS *interface_faces_is)
1262*c4762a1bSJed Brown {
1263*c4762a1bSJed Brown   MPI_Comm               comm;
1264*c4762a1bSJed Brown   DMLabel                label;
1265*c4762a1bSJed Brown   IS                     part_boundary_faces_is;
1266*c4762a1bSJed Brown   const char             partBoundaryName[] = "DMPlexDistributeInterpolateMarkInterface_partBoundary";
1267*c4762a1bSJed Brown   PetscInt               value = 1;
1268*c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1269*c4762a1bSJed Brown   PetscErrorCode         ierr;
1270*c4762a1bSJed Brown 
1271*c4762a1bSJed Brown   PetscFunctionBegin;
1272*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)ipdm, &comm);CHKERRQ(ierr);
1273*c4762a1bSJed Brown   ierr = DMPlexIsInterpolatedCollective(ipdm, &intp);CHKERRQ(ierr);
1274*c4762a1bSJed Brown   if (intp != DMPLEX_INTERPOLATED_FULL) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex");
1275*c4762a1bSJed Brown 
1276*c4762a1bSJed Brown   /* get ipdm partition boundary (partBoundary) */
1277*c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, partBoundaryName, &label);CHKERRQ(ierr);
1278*c4762a1bSJed Brown   ierr = DMAddLabel(ipdm, label);CHKERRQ(ierr);
1279*c4762a1bSJed Brown   ierr = DMPlexMarkBoundaryFaces(ipdm, value, label);CHKERRQ(ierr);
1280*c4762a1bSJed Brown   ierr = DMLabelViewFromOptionsOnComm_Private(label, "-ipdm_part_boundary_view", comm);CHKERRQ(ierr);
1281*c4762a1bSJed Brown   ierr = DMLabelGetStratumISOnComm_Private(label, value, comm, &part_boundary_faces_is);CHKERRQ(ierr);
1282*c4762a1bSJed Brown   ierr = DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE);CHKERRQ(ierr);
1283*c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
1284*c4762a1bSJed Brown 
1285*c4762a1bSJed Brown   /* remove ipdm whole-domain boundary (boundary_faces_is) from ipdm partition boundary (part_boundary_faces_is), resulting just in inter-partition interface */
1286*c4762a1bSJed Brown   ierr = ISDifference(part_boundary_faces_is,boundary_faces_is,interface_faces_is);CHKERRQ(ierr);
1287*c4762a1bSJed Brown   ierr = ISDestroy(&part_boundary_faces_is);CHKERRQ(ierr);
1288*c4762a1bSJed Brown   PetscFunctionReturn(0);
1289*c4762a1bSJed Brown }
1290*c4762a1bSJed Brown 
1291*c4762a1bSJed Brown /* compute inter-partition interface including edges and vertices */
1292*c4762a1bSJed Brown static PetscErrorCode DMPlexComputeCompleteInterface_Private(DM ipdm, IS interface_faces_is, IS *interface_is)
1293*c4762a1bSJed Brown {
1294*c4762a1bSJed Brown   DMLabel                label;
1295*c4762a1bSJed Brown   PetscInt               value = 1;
1296*c4762a1bSJed Brown   const char             interfaceName[] = "DMPlexDistributeInterpolateMarkInterface_interface";
1297*c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1298*c4762a1bSJed Brown   MPI_Comm               comm;
1299*c4762a1bSJed Brown   PetscErrorCode         ierr;
1300*c4762a1bSJed Brown 
1301*c4762a1bSJed Brown   PetscFunctionBegin;
1302*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)ipdm, &comm);CHKERRQ(ierr);
1303*c4762a1bSJed Brown   ierr = DMPlexIsInterpolatedCollective(ipdm, &intp);CHKERRQ(ierr);
1304*c4762a1bSJed Brown   if (intp != DMPLEX_INTERPOLATED_FULL) SETERRQ(comm, PETSC_ERR_ARG_WRONG, "only for fully interpolated DMPlex");
1305*c4762a1bSJed Brown 
1306*c4762a1bSJed Brown   ierr = DMLabelCreate(PETSC_COMM_SELF, interfaceName, &label);CHKERRQ(ierr);
1307*c4762a1bSJed Brown   ierr = DMAddLabel(ipdm, label);CHKERRQ(ierr);
1308*c4762a1bSJed Brown   ierr = DMLabelSetStratumIS(label, value, interface_faces_is);CHKERRQ(ierr);
1309*c4762a1bSJed Brown   ierr = DMLabelViewFromOptionsOnComm_Private(label, "-interface_faces_view", comm);CHKERRQ(ierr);
1310*c4762a1bSJed Brown   ierr = DMPlexLabelComplete(ipdm, label);CHKERRQ(ierr);
1311*c4762a1bSJed Brown   ierr = DMLabelViewFromOptionsOnComm_Private(label, "-interface_view", comm);CHKERRQ(ierr);
1312*c4762a1bSJed Brown   ierr = DMLabelGetStratumISOnComm_Private(label, value, comm, interface_is);CHKERRQ(ierr);
1313*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)*interface_is, "interface_is");CHKERRQ(ierr);
1314*c4762a1bSJed Brown   ierr = ISViewFromOptions(*interface_is, NULL, "-interface_is_view");CHKERRQ(ierr);
1315*c4762a1bSJed Brown   ierr = DMRemoveLabelBySelf(ipdm, &label, PETSC_TRUE);CHKERRQ(ierr);
1316*c4762a1bSJed Brown   ierr = DMLabelDestroy(&label);CHKERRQ(ierr);
1317*c4762a1bSJed Brown   PetscFunctionReturn(0);
1318*c4762a1bSJed Brown }
1319*c4762a1bSJed Brown 
1320*c4762a1bSJed Brown static PetscErrorCode PointSFGetOutwardInterfacePoints(PetscSF sf, IS *is)
1321*c4762a1bSJed Brown {
1322*c4762a1bSJed Brown   PetscInt        n;
1323*c4762a1bSJed Brown   const PetscInt  *arr;
1324*c4762a1bSJed Brown   PetscErrorCode  ierr;
1325*c4762a1bSJed Brown 
1326*c4762a1bSJed Brown   PetscFunctionBegin;
1327*c4762a1bSJed Brown   ierr = PetscSFGetGraph(sf, NULL, &n, &arr, NULL);CHKERRQ(ierr);
1328*c4762a1bSJed Brown   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_USE_POINTER, is);CHKERRQ(ierr);
1329*c4762a1bSJed Brown   PetscFunctionReturn(0);
1330*c4762a1bSJed Brown }
1331*c4762a1bSJed Brown 
1332*c4762a1bSJed Brown static PetscErrorCode PointSFGetInwardInterfacePoints(PetscSF sf, IS *is)
1333*c4762a1bSJed Brown {
1334*c4762a1bSJed Brown   PetscInt        n;
1335*c4762a1bSJed Brown   const PetscInt  *rootdegree;
1336*c4762a1bSJed Brown   PetscInt        *arr;
1337*c4762a1bSJed Brown   PetscErrorCode  ierr;
1338*c4762a1bSJed Brown 
1339*c4762a1bSJed Brown   PetscFunctionBegin;
1340*c4762a1bSJed Brown   ierr = PetscSFSetUp(sf);CHKERRQ(ierr);
1341*c4762a1bSJed Brown   ierr = PetscSFComputeDegreeBegin(sf, &rootdegree);CHKERRQ(ierr);
1342*c4762a1bSJed Brown   ierr = PetscSFComputeDegreeEnd(sf, &rootdegree);CHKERRQ(ierr);
1343*c4762a1bSJed Brown   ierr = PetscSFComputeMultiRootOriginalNumbering(sf, rootdegree, &n, &arr);CHKERRQ(ierr);
1344*c4762a1bSJed Brown   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)sf), n, arr, PETSC_OWN_POINTER, is);CHKERRQ(ierr);
1345*c4762a1bSJed Brown   PetscFunctionReturn(0);
1346*c4762a1bSJed Brown }
1347*c4762a1bSJed Brown 
1348*c4762a1bSJed Brown static PetscErrorCode PointSFGetInterfacePoints_Private(PetscSF pointSF, IS *is)
1349*c4762a1bSJed Brown {
1350*c4762a1bSJed Brown   IS pointSF_out_is, pointSF_in_is;
1351*c4762a1bSJed Brown   PetscErrorCode  ierr;
1352*c4762a1bSJed Brown 
1353*c4762a1bSJed Brown   PetscFunctionBegin;
1354*c4762a1bSJed Brown   ierr = PointSFGetOutwardInterfacePoints(pointSF, &pointSF_out_is);CHKERRQ(ierr);
1355*c4762a1bSJed Brown   ierr = PointSFGetInwardInterfacePoints(pointSF, &pointSF_in_is);CHKERRQ(ierr);
1356*c4762a1bSJed Brown   ierr = ISExpand(pointSF_out_is, pointSF_in_is, is);CHKERRQ(ierr);
1357*c4762a1bSJed Brown   ierr = ISDestroy(&pointSF_out_is);CHKERRQ(ierr);
1358*c4762a1bSJed Brown   ierr = ISDestroy(&pointSF_in_is);CHKERRQ(ierr);
1359*c4762a1bSJed Brown   PetscFunctionReturn(0);
1360*c4762a1bSJed Brown }
1361*c4762a1bSJed Brown 
1362*c4762a1bSJed Brown #define MYCHKERRQ(ierr) do {if (PetscUnlikely(ierr)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "PointSF is wrong. Unable to show details!");} while (0)
1363*c4762a1bSJed Brown 
1364*c4762a1bSJed Brown static PetscErrorCode ViewPointsWithType_Internal(DM dm, IS pointsIS, PetscViewer v)
1365*c4762a1bSJed Brown {
1366*c4762a1bSJed Brown   DMLabel         label;
1367*c4762a1bSJed Brown   PetscSection    coordsSection;
1368*c4762a1bSJed Brown   Vec             coordsVec;
1369*c4762a1bSJed Brown   PetscScalar     *coordsScalar;
1370*c4762a1bSJed Brown   PetscInt        coneSize, depth, dim, i, p, npoints;
1371*c4762a1bSJed Brown   const PetscInt  *points;
1372*c4762a1bSJed Brown   PetscErrorCode  ierr;
1373*c4762a1bSJed Brown 
1374*c4762a1bSJed Brown   PetscFunctionBegin;
1375*c4762a1bSJed Brown   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1376*c4762a1bSJed Brown   ierr = DMGetCoordinateSection(dm, &coordsSection);CHKERRQ(ierr);
1377*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(dm, &coordsVec);CHKERRQ(ierr);
1378*c4762a1bSJed Brown   ierr = VecGetArray(coordsVec, &coordsScalar);CHKERRQ(ierr);
1379*c4762a1bSJed Brown   ierr = ISGetLocalSize(pointsIS, &npoints);CHKERRQ(ierr);
1380*c4762a1bSJed Brown   ierr = ISGetIndices(pointsIS, &points);CHKERRQ(ierr);
1381*c4762a1bSJed Brown   ierr = DMPlexGetDepthLabel(dm, &label);CHKERRQ(ierr);
1382*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr);
1383*c4762a1bSJed Brown   for (i=0; i<npoints; i++) {
1384*c4762a1bSJed Brown     p = points[i];
1385*c4762a1bSJed Brown     ierr = DMLabelGetValue(label, p, &depth);CHKERRQ(ierr);
1386*c4762a1bSJed Brown     if (!depth) {
1387*c4762a1bSJed Brown       PetscInt        c, n, o;
1388*c4762a1bSJed Brown       PetscReal       coords[3];
1389*c4762a1bSJed Brown       char            coordstr[128];
1390*c4762a1bSJed Brown 
1391*c4762a1bSJed Brown       ierr = PetscSectionGetDof(coordsSection, p, &n);CHKERRQ(ierr);
1392*c4762a1bSJed Brown       ierr = PetscSectionGetOffset(coordsSection, p, &o);CHKERRQ(ierr);
1393*c4762a1bSJed Brown       for (c=0; c<n; c++) coords[c] = PetscRealPart(coordsScalar[o+c]);
1394*c4762a1bSJed Brown       ierr = coord2str(coordstr, sizeof(coordstr), n, coords, 1.0);CHKERRQ(ierr);
1395*c4762a1bSJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(v, "vertex %D w/ coordinates %s\n", p, coordstr);CHKERRQ(ierr);
1396*c4762a1bSJed Brown     } else {
1397*c4762a1bSJed Brown       char            entityType[16];
1398*c4762a1bSJed Brown 
1399*c4762a1bSJed Brown       switch (depth) {
1400*c4762a1bSJed Brown         case 1: ierr = PetscStrcpy(entityType, "edge");CHKERRQ(ierr); break;
1401*c4762a1bSJed Brown         case 2: ierr = PetscStrcpy(entityType, "face");CHKERRQ(ierr); break;
1402*c4762a1bSJed Brown         case 3: ierr = PetscStrcpy(entityType, "cell");CHKERRQ(ierr); break;
1403*c4762a1bSJed Brown         default: SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for depth <= 3");CHKERRQ(ierr);
1404*c4762a1bSJed Brown       }
1405*c4762a1bSJed Brown       if (depth == dim && dim < 3) {
1406*c4762a1bSJed Brown         ierr = PetscStrlcat(entityType, " (cell)", sizeof(entityType));CHKERRQ(ierr);
1407*c4762a1bSJed Brown       }
1408*c4762a1bSJed Brown       ierr = PetscViewerASCIISynchronizedPrintf(v, "%s %D\n", entityType, p);CHKERRQ(ierr);
1409*c4762a1bSJed Brown     }
1410*c4762a1bSJed Brown     ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
1411*c4762a1bSJed Brown     if (coneSize) {
1412*c4762a1bSJed Brown       const PetscInt *cone;
1413*c4762a1bSJed Brown       IS             coneIS;
1414*c4762a1bSJed Brown 
1415*c4762a1bSJed Brown       ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1416*c4762a1bSJed Brown       ierr = ISCreateGeneral(PETSC_COMM_SELF, coneSize, cone, PETSC_USE_POINTER, &coneIS);CHKERRQ(ierr);
1417*c4762a1bSJed Brown       ierr = ViewPointsWithType_Internal(dm, coneIS, v);CHKERRQ(ierr);
1418*c4762a1bSJed Brown       ierr = ISDestroy(&coneIS);CHKERRQ(ierr);
1419*c4762a1bSJed Brown     }
1420*c4762a1bSJed Brown   }
1421*c4762a1bSJed Brown   ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr);
1422*c4762a1bSJed Brown   ierr = VecRestoreArray(coordsVec, &coordsScalar);CHKERRQ(ierr);
1423*c4762a1bSJed Brown   ierr = ISRestoreIndices(pointsIS, &points);CHKERRQ(ierr);
1424*c4762a1bSJed Brown   PetscFunctionReturn(0);
1425*c4762a1bSJed Brown }
1426*c4762a1bSJed Brown 
1427*c4762a1bSJed Brown static PetscErrorCode ViewPointsWithType(DM dm, IS points, PetscViewer v)
1428*c4762a1bSJed Brown {
1429*c4762a1bSJed Brown   PetscBool       flg;
1430*c4762a1bSJed Brown   PetscInt        npoints;
1431*c4762a1bSJed Brown   PetscMPIInt     rank;
1432*c4762a1bSJed Brown   PetscErrorCode  ierr;
1433*c4762a1bSJed Brown 
1434*c4762a1bSJed Brown   PetscFunctionBegin;
1435*c4762a1bSJed Brown   ierr = PetscObjectTypeCompare((PetscObject)v, PETSCVIEWERASCII, &flg);CHKERRQ(ierr);
1436*c4762a1bSJed Brown   if (!flg) SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Only for ASCII viewer");CHKERRQ(ierr);
1437*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)v), &rank);CHKERRQ(ierr);
1438*c4762a1bSJed Brown   ierr = PetscViewerASCIIPushSynchronized(v);CHKERRQ(ierr);
1439*c4762a1bSJed Brown   ierr = ISGetLocalSize(points, &npoints);CHKERRQ(ierr);
1440*c4762a1bSJed Brown   if (npoints) {
1441*c4762a1bSJed Brown     ierr = PetscViewerASCIISynchronizedPrintf(v, "[%d] --------\n", rank);CHKERRQ(ierr);
1442*c4762a1bSJed Brown     ierr = ViewPointsWithType_Internal(dm, points, v);CHKERRQ(ierr);
1443*c4762a1bSJed Brown   }
1444*c4762a1bSJed Brown   ierr = PetscViewerFlush(v);CHKERRQ(ierr);
1445*c4762a1bSJed Brown   ierr = PetscViewerASCIIPopSynchronized(v);CHKERRQ(ierr);
1446*c4762a1bSJed Brown   PetscFunctionReturn(0);
1447*c4762a1bSJed Brown }
1448*c4762a1bSJed Brown 
1449*c4762a1bSJed Brown static PetscErrorCode DMPlexComparePointSFWithInterface_Private(DM ipdm, IS interface_is)
1450*c4762a1bSJed Brown {
1451*c4762a1bSJed Brown   PetscSF         pointsf;
1452*c4762a1bSJed Brown   IS              pointsf_is;
1453*c4762a1bSJed Brown   PetscBool       flg;
1454*c4762a1bSJed Brown   MPI_Comm        comm;
1455*c4762a1bSJed Brown   PetscMPIInt     size;
1456*c4762a1bSJed Brown   PetscErrorCode  ierr;
1457*c4762a1bSJed Brown 
1458*c4762a1bSJed Brown   PetscFunctionBegin;
1459*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)ipdm, &comm);CHKERRQ(ierr);
1460*c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1461*c4762a1bSJed Brown   ierr = DMGetPointSF(ipdm, &pointsf);CHKERRQ(ierr);
1462*c4762a1bSJed Brown   if (pointsf) {
1463*c4762a1bSJed Brown     PetscInt nroots;
1464*c4762a1bSJed Brown     ierr = PetscSFGetGraph(pointsf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1465*c4762a1bSJed Brown     if (nroots < 0) pointsf = NULL; /* uninitialized SF */
1466*c4762a1bSJed Brown   }
1467*c4762a1bSJed Brown   if (!pointsf) {
1468*c4762a1bSJed Brown     PetscInt N=0;
1469*c4762a1bSJed Brown     if (interface_is) {ierr = ISGetSize(interface_is, &N);CHKERRQ(ierr);}
1470*c4762a1bSJed Brown     if (N) SETERRQ(comm, PETSC_ERR_PLIB, "interface_is should be NULL or empty for PointSF being NULL");
1471*c4762a1bSJed Brown     PetscFunctionReturn(0);
1472*c4762a1bSJed Brown   }
1473*c4762a1bSJed Brown 
1474*c4762a1bSJed Brown   /* get PointSF points as IS pointsf_is */
1475*c4762a1bSJed Brown   ierr = PointSFGetInterfacePoints_Private(pointsf, &pointsf_is);CHKERRQ(ierr);
1476*c4762a1bSJed Brown 
1477*c4762a1bSJed Brown   /* compare pointsf_is with interface_is */
1478*c4762a1bSJed Brown   ierr = ISEqual(interface_is, pointsf_is, &flg);CHKERRQ(ierr);
1479*c4762a1bSJed Brown   ierr = MPI_Allreduce(MPI_IN_PLACE,&flg,1,MPIU_BOOL,MPI_LAND,comm);CHKERRQ(ierr);
1480*c4762a1bSJed Brown   if (!flg) {
1481*c4762a1bSJed Brown     IS pointsf_extra_is, pointsf_missing_is;
1482*c4762a1bSJed Brown     PetscViewer errv = PETSC_VIEWER_STDERR_(comm);
1483*c4762a1bSJed Brown     ierr = ISDifference(interface_is, pointsf_is, &pointsf_missing_is);MYCHKERRQ(ierr);
1484*c4762a1bSJed Brown     ierr = ISDifference(pointsf_is, interface_is, &pointsf_extra_is);MYCHKERRQ(ierr);
1485*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(errv, "Points missing in PointSF:\n");MYCHKERRQ(ierr);
1486*c4762a1bSJed Brown     ierr = ViewPointsWithType(ipdm, pointsf_missing_is, errv);MYCHKERRQ(ierr);
1487*c4762a1bSJed Brown     ierr = PetscViewerASCIIPrintf(errv, "Extra points in PointSF:\n");MYCHKERRQ(ierr);
1488*c4762a1bSJed Brown     ierr = ViewPointsWithType(ipdm, pointsf_extra_is, errv);MYCHKERRQ(ierr);
1489*c4762a1bSJed Brown     ierr = ISDestroy(&pointsf_extra_is);MYCHKERRQ(ierr);
1490*c4762a1bSJed Brown     ierr = ISDestroy(&pointsf_missing_is);MYCHKERRQ(ierr);
1491*c4762a1bSJed Brown     SETERRQ(comm, PETSC_ERR_PLIB, "PointSF is wrong! See details above.");
1492*c4762a1bSJed Brown   }
1493*c4762a1bSJed Brown   ierr = ISDestroy(&pointsf_is);CHKERRQ(ierr);
1494*c4762a1bSJed Brown   PetscFunctionReturn(0);
1495*c4762a1bSJed Brown }
1496*c4762a1bSJed Brown 
1497*c4762a1bSJed Brown /* remove faces & edges from label, leave just vertices */
1498*c4762a1bSJed Brown static PetscErrorCode DMPlexISFilterVertices_Private(DM dm, IS points)
1499*c4762a1bSJed Brown {
1500*c4762a1bSJed Brown   PetscInt        vStart, vEnd;
1501*c4762a1bSJed Brown   MPI_Comm        comm;
1502*c4762a1bSJed Brown   PetscErrorCode  ierr;
1503*c4762a1bSJed Brown 
1504*c4762a1bSJed Brown   PetscFunctionBegin;
1505*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1506*c4762a1bSJed Brown   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
1507*c4762a1bSJed Brown   ierr = ISGeneralFilter(points, vStart, vEnd);CHKERRQ(ierr);
1508*c4762a1bSJed Brown   PetscFunctionReturn(0);
1509*c4762a1bSJed Brown }
1510*c4762a1bSJed Brown 
1511*c4762a1bSJed Brown /*
1512*c4762a1bSJed Brown   DMPlexCheckPointSFHeavy - Thorough test that the PointSF after parallel DMPlexInterpolate() includes exactly all interface points.
1513*c4762a1bSJed Brown 
1514*c4762a1bSJed Brown   Collective
1515*c4762a1bSJed Brown 
1516*c4762a1bSJed Brown   Input Parameters:
1517*c4762a1bSJed Brown . dm - The DMPlex object
1518*c4762a1bSJed Brown 
1519*c4762a1bSJed Brown   Notes:
1520*c4762a1bSJed Brown   The input DMPlex must be serial (one partition has all points, the other partitions have no points).
1521*c4762a1bSJed 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).
1522*c4762a1bSJed Brown   This is mainly intended for debugging/testing purposes.
1523*c4762a1bSJed Brown 
1524*c4762a1bSJed Brown   Algorithm:
1525*c4762a1bSJed Brown   1. boundary faces of the serial version of the whole mesh are found using DMPlexMarkBoundaryFaces()
1526*c4762a1bSJed 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
1527*c4762a1bSJed Brown   3. the mesh is distributed or loaded in parallel
1528*c4762a1bSJed Brown   4. boundary faces of the distributed mesh are reconstructed from PortableBoundary using DMPlexFindVertices()
1529*c4762a1bSJed Brown   5. partition boundary faces of the parallel mesh are found using DMPlexMarkBoundaryFaces()
1530*c4762a1bSJed Brown   6. partition interfaces are computed as set difference of partition boundary faces minus the reconstructed boundary
1531*c4762a1bSJed 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
1532*c4762a1bSJed Brown 
1533*c4762a1bSJed Brown   Level: developer
1534*c4762a1bSJed Brown 
1535*c4762a1bSJed Brown .seealso: DMGetPointSF(), DMPlexCheckSymmetry(), DMPlexCheckSkeleton(), DMPlexCheckFaces()
1536*c4762a1bSJed Brown */
1537*c4762a1bSJed Brown static PetscErrorCode DMPlexCheckPointSFHeavy(DM dm, PortableBoundary bnd)
1538*c4762a1bSJed Brown {
1539*c4762a1bSJed Brown   DM                     ipdm=NULL;
1540*c4762a1bSJed Brown   IS                     boundary_faces_is, interface_faces_is, interface_is;
1541*c4762a1bSJed Brown   DMPlexInterpolatedFlag intp;
1542*c4762a1bSJed Brown   MPI_Comm               comm;
1543*c4762a1bSJed Brown   PetscErrorCode         ierr;
1544*c4762a1bSJed Brown 
1545*c4762a1bSJed Brown   PetscFunctionBegin;
1546*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1547*c4762a1bSJed Brown 
1548*c4762a1bSJed Brown   ierr = DMPlexIsInterpolatedCollective(dm, &intp);CHKERRQ(ierr);
1549*c4762a1bSJed Brown   if (intp == DMPLEX_INTERPOLATED_FULL) {
1550*c4762a1bSJed Brown     ipdm = dm;
1551*c4762a1bSJed Brown   } else {
1552*c4762a1bSJed Brown     /* create temporary interpolated DM if input DM is not interpolated */
1553*c4762a1bSJed Brown     ierr = DMPlexSetOrientInterface_Private(dm, PETSC_FALSE);CHKERRQ(ierr);
1554*c4762a1bSJed Brown     ierr = DMPlexInterpolate(dm, &ipdm);CHKERRQ(ierr); /* with DMPlexOrientInterface_Internal() call skipped so that PointSF issues are left to DMPlexComparePointSFWithInterface_Private() below */
1555*c4762a1bSJed Brown     ierr = DMPlexSetOrientInterface_Private(dm, PETSC_TRUE);CHKERRQ(ierr);
1556*c4762a1bSJed Brown   }
1557*c4762a1bSJed Brown   ierr = DMViewFromOptions(ipdm, NULL, "-ipdm_view");CHKERRQ(ierr);
1558*c4762a1bSJed Brown 
1559*c4762a1bSJed Brown   /* recover ipdm whole-domain boundary faces from the expanded vertices coordinates */
1560*c4762a1bSJed Brown   ierr = DMPlexExpandedVerticesCoordinatesToFaces_Private(ipdm, bnd, &boundary_faces_is);CHKERRQ(ierr);
1561*c4762a1bSJed Brown   /* get inter-partition interface faces (interface_faces_is)*/
1562*c4762a1bSJed Brown   ierr = DMPlexGetInterfaceFaces_Private(ipdm, boundary_faces_is, &interface_faces_is);CHKERRQ(ierr);
1563*c4762a1bSJed Brown   /* compute inter-partition interface including edges and vertices (interface_is) */
1564*c4762a1bSJed Brown   ierr = DMPlexComputeCompleteInterface_Private(ipdm, interface_faces_is, &interface_is);CHKERRQ(ierr);
1565*c4762a1bSJed Brown   /* destroy immediate ISs */
1566*c4762a1bSJed Brown   ierr = ISDestroy(&boundary_faces_is);CHKERRQ(ierr);
1567*c4762a1bSJed Brown   ierr = ISDestroy(&interface_faces_is);CHKERRQ(ierr);
1568*c4762a1bSJed Brown 
1569*c4762a1bSJed Brown   /* for uninterpolated case, keep just vertices in interface */
1570*c4762a1bSJed Brown   if (!intp) {
1571*c4762a1bSJed Brown     ierr = DMPlexISFilterVertices_Private(ipdm, interface_is);CHKERRQ(ierr);
1572*c4762a1bSJed Brown     ierr = DMDestroy(&ipdm);CHKERRQ(ierr);
1573*c4762a1bSJed Brown   }
1574*c4762a1bSJed Brown 
1575*c4762a1bSJed Brown   /* compare PointSF with the boundary reconstructed from coordinates */
1576*c4762a1bSJed Brown   ierr = DMPlexComparePointSFWithInterface_Private(dm, interface_is);CHKERRQ(ierr);
1577*c4762a1bSJed Brown   ierr = PetscPrintf(comm, "DMPlexCheckPointSFHeavy PASSED\n");CHKERRQ(ierr);
1578*c4762a1bSJed Brown   ierr = ISDestroy(&interface_is);CHKERRQ(ierr);
1579*c4762a1bSJed Brown   PetscFunctionReturn(0);
1580*c4762a1bSJed Brown }
1581*c4762a1bSJed Brown 
1582*c4762a1bSJed Brown int main(int argc, char **argv)
1583*c4762a1bSJed Brown {
1584*c4762a1bSJed Brown   DM             dm;
1585*c4762a1bSJed Brown   AppCtx         user;
1586*c4762a1bSJed Brown   PetscErrorCode ierr;
1587*c4762a1bSJed Brown 
1588*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
1589*c4762a1bSJed Brown   ierr = PetscLogStageRegister("create",&stage[0]);CHKERRQ(ierr);
1590*c4762a1bSJed Brown   ierr = PetscLogStageRegister("distribute",&stage[1]);CHKERRQ(ierr);
1591*c4762a1bSJed Brown   ierr = PetscLogStageRegister("interpolate",&stage[2]);CHKERRQ(ierr);
1592*c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
1593*c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
1594*c4762a1bSJed Brown   if (user.nPointsToExpand) {
1595*c4762a1bSJed Brown     ierr = TestExpandPoints(dm, &user);CHKERRQ(ierr);
1596*c4762a1bSJed Brown   }
1597*c4762a1bSJed Brown   if (user.ncoords) {
1598*c4762a1bSJed Brown     ierr = ViewVerticesFromCoords(dm, user.ncoords/user.dim, user.coords, user.coordsTol, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1599*c4762a1bSJed Brown   }
1600*c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
1601*c4762a1bSJed Brown   ierr = PetscFinalize();
1602*c4762a1bSJed Brown   return ierr;
1603*c4762a1bSJed Brown }
1604*c4762a1bSJed Brown 
1605*c4762a1bSJed Brown /*TEST
1606*c4762a1bSJed Brown 
1607*c4762a1bSJed Brown   testset:
1608*c4762a1bSJed Brown     nsize: 2
1609*c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail
1610*c4762a1bSJed Brown     args: -dm_plex_check_all
1611*c4762a1bSJed Brown     test:
1612*c4762a1bSJed Brown       suffix: 1_tri_dist0
1613*c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1614*c4762a1bSJed Brown     test:
1615*c4762a1bSJed Brown       suffix: 1_tri_dist1
1616*c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1617*c4762a1bSJed Brown     test:
1618*c4762a1bSJed Brown       suffix: 1_quad_dist0
1619*c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output}
1620*c4762a1bSJed Brown     test:
1621*c4762a1bSJed Brown       suffix: 1_quad_dist1
1622*c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output}
1623*c4762a1bSJed Brown     test:
1624*c4762a1bSJed Brown       suffix: 1_1d_dist1
1625*c4762a1bSJed Brown       args: -dim 1 -distribute 1
1626*c4762a1bSJed Brown 
1627*c4762a1bSJed Brown   testset:
1628*c4762a1bSJed Brown     nsize: 3
1629*c4762a1bSJed Brown     args: -testnum 1 -interpolate create
1630*c4762a1bSJed Brown     args: -dm_plex_check_all
1631*c4762a1bSJed Brown     test:
1632*c4762a1bSJed Brown       suffix: 2
1633*c4762a1bSJed Brown       args: -dm_view ascii::ascii_info_detail
1634*c4762a1bSJed Brown     test:
1635*c4762a1bSJed Brown       suffix: 2a
1636*c4762a1bSJed Brown       args: -dm_plex_check_cones_conform_on_interfaces_verbose
1637*c4762a1bSJed Brown     test:
1638*c4762a1bSJed Brown       suffix: 2b
1639*c4762a1bSJed Brown       args: -test_expand_points 0,1,2,5,6
1640*c4762a1bSJed Brown     test:
1641*c4762a1bSJed Brown       suffix: 2c
1642*c4762a1bSJed Brown       args: -test_expand_points 0,1,2,5,6 -test_expand_points_empty
1643*c4762a1bSJed Brown 
1644*c4762a1bSJed Brown   testset:
1645*c4762a1bSJed Brown     # the same as 1% for 3D
1646*c4762a1bSJed Brown     nsize: 2
1647*c4762a1bSJed Brown     args: -dim 3 -dm_view ascii::ascii_info_detail
1648*c4762a1bSJed Brown     args: -dm_plex_check_all
1649*c4762a1bSJed Brown     test:
1650*c4762a1bSJed Brown       suffix: 4_tet_dist0
1651*c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1652*c4762a1bSJed Brown     test:
1653*c4762a1bSJed Brown       suffix: 4_tet_dist1
1654*c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1655*c4762a1bSJed Brown     test:
1656*c4762a1bSJed Brown       suffix: 4_hex_dist0
1657*c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 0 -interpolate {{none create}separate output}
1658*c4762a1bSJed Brown     test:
1659*c4762a1bSJed Brown       suffix: 4_hex_dist1
1660*c4762a1bSJed Brown       args: -cell_simplex 0 -distribute 1 -interpolate {{none create after_distribute}separate output}
1661*c4762a1bSJed Brown 
1662*c4762a1bSJed Brown   test:
1663*c4762a1bSJed Brown     # the same as 4_tet_dist0 but test different initial orientations
1664*c4762a1bSJed Brown     suffix: 4_tet_test_orient
1665*c4762a1bSJed Brown     nsize: 2
1666*c4762a1bSJed Brown     args: -dim 3 -distribute 0
1667*c4762a1bSJed Brown     args: -dm_plex_check_all
1668*c4762a1bSJed Brown     args: -rotate_interface_0 {{0 1 2 11 12 13}}
1669*c4762a1bSJed Brown     args: -rotate_interface_1 {{0 1 2 11 12 13}}
1670*c4762a1bSJed Brown 
1671*c4762a1bSJed Brown   testset:
1672*c4762a1bSJed Brown     requires: exodusii
1673*c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/TwoQuads.exo
1674*c4762a1bSJed Brown     args: -dm_view ascii::ascii_info_detail
1675*c4762a1bSJed Brown     args: -dm_plex_check_all
1676*c4762a1bSJed Brown     args: -custom_view
1677*c4762a1bSJed Brown     test:
1678*c4762a1bSJed Brown       suffix: 5_seq
1679*c4762a1bSJed Brown       nsize: 1
1680*c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1681*c4762a1bSJed Brown     test:
1682*c4762a1bSJed Brown       suffix: 5_dist0
1683*c4762a1bSJed Brown       nsize: 2
1684*c4762a1bSJed Brown       args: -distribute 0 -interpolate {{none create}separate output}
1685*c4762a1bSJed Brown     test:
1686*c4762a1bSJed Brown       suffix: 5_dist1
1687*c4762a1bSJed Brown       nsize: 2
1688*c4762a1bSJed Brown       args: -distribute 1 -interpolate {{none create after_distribute}separate output}
1689*c4762a1bSJed Brown 
1690*c4762a1bSJed Brown   testset:
1691*c4762a1bSJed Brown     nsize: {{1 2 4}}
1692*c4762a1bSJed Brown     args: -use_generator
1693*c4762a1bSJed Brown     args: -dm_plex_check_all
1694*c4762a1bSJed Brown     args: -distribute -interpolate none
1695*c4762a1bSJed Brown     test:
1696*c4762a1bSJed Brown       suffix: 6_tri
1697*c4762a1bSJed Brown       requires: triangle
1698*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_plex_generator triangle
1699*c4762a1bSJed Brown     test:
1700*c4762a1bSJed Brown       suffix: 6_quad
1701*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1702*c4762a1bSJed Brown     test:
1703*c4762a1bSJed Brown       suffix: 6_tet
1704*c4762a1bSJed Brown       requires: ctetgen
1705*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_plex_generator ctetgen
1706*c4762a1bSJed Brown     test:
1707*c4762a1bSJed Brown       suffix: 6_hex
1708*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1709*c4762a1bSJed Brown   testset:
1710*c4762a1bSJed Brown     nsize: {{1 2 4}}
1711*c4762a1bSJed Brown     args: -use_generator
1712*c4762a1bSJed Brown     args: -dm_plex_check_all
1713*c4762a1bSJed Brown     args: -distribute -interpolate create
1714*c4762a1bSJed Brown     test:
1715*c4762a1bSJed Brown       suffix: 6_int_tri
1716*c4762a1bSJed Brown       requires: triangle
1717*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_plex_generator triangle
1718*c4762a1bSJed Brown     test:
1719*c4762a1bSJed Brown       suffix: 6_int_quad
1720*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1721*c4762a1bSJed Brown     test:
1722*c4762a1bSJed Brown       suffix: 6_int_tet
1723*c4762a1bSJed Brown       requires: ctetgen
1724*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_plex_generator ctetgen
1725*c4762a1bSJed Brown     test:
1726*c4762a1bSJed Brown       suffix: 6_int_hex
1727*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1728*c4762a1bSJed Brown   testset:
1729*c4762a1bSJed Brown     nsize: {{2 4}}
1730*c4762a1bSJed Brown     args: -use_generator
1731*c4762a1bSJed Brown     args: -dm_plex_check_all
1732*c4762a1bSJed Brown     args: -distribute -interpolate after_distribute
1733*c4762a1bSJed Brown     test:
1734*c4762a1bSJed Brown       suffix: 6_parint_tri
1735*c4762a1bSJed Brown       requires: triangle
1736*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 1 -dm_plex_generator triangle
1737*c4762a1bSJed Brown     test:
1738*c4762a1bSJed Brown       suffix: 6_parint_quad
1739*c4762a1bSJed Brown       args: -faces {{2,2  1,3  7,4}} -cell_simplex 0
1740*c4762a1bSJed Brown     test:
1741*c4762a1bSJed Brown       suffix: 6_parint_tet
1742*c4762a1bSJed Brown       requires: ctetgen
1743*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 1 -dm_plex_generator ctetgen
1744*c4762a1bSJed Brown     test:
1745*c4762a1bSJed Brown       suffix: 6_parint_hex
1746*c4762a1bSJed Brown       args: -faces {{2,2,2  1,3,5  3,4,7}} -cell_simplex 0
1747*c4762a1bSJed Brown 
1748*c4762a1bSJed Brown   testset: # 7 EXODUS
1749*c4762a1bSJed Brown     requires: exodusii
1750*c4762a1bSJed Brown     args: -dm_plex_check_all
1751*c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo
1752*c4762a1bSJed Brown     args: -distribute
1753*c4762a1bSJed Brown     test: # seq load, simple partitioner
1754*c4762a1bSJed Brown       suffix: 7_exo
1755*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1756*c4762a1bSJed Brown       args: -interpolate none
1757*c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1758*c4762a1bSJed Brown       suffix: 7_exo_int_simple
1759*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1760*c4762a1bSJed Brown       args: -interpolate create
1761*c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1762*c4762a1bSJed Brown       suffix: 7_exo_int_metis
1763*c4762a1bSJed Brown       requires: parmetis
1764*c4762a1bSJed Brown       nsize: {{2 4 5}}
1765*c4762a1bSJed Brown       args: -interpolate create
1766*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1767*c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1768*c4762a1bSJed Brown       suffix: 7_exo_simple_int
1769*c4762a1bSJed Brown       nsize: {{2 4 5}}
1770*c4762a1bSJed Brown       args: -interpolate after_distribute
1771*c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1772*c4762a1bSJed Brown       suffix: 7_exo_metis_int
1773*c4762a1bSJed Brown       requires: parmetis
1774*c4762a1bSJed Brown       nsize: {{2 4 5}}
1775*c4762a1bSJed Brown       args: -interpolate after_distribute
1776*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1777*c4762a1bSJed Brown 
1778*c4762a1bSJed Brown   testset: # 7 HDF5 SEQUANTIAL LOAD
1779*c4762a1bSJed Brown     requires: hdf5 !complex
1780*c4762a1bSJed Brown     args: -dm_plex_check_all
1781*c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1782*c4762a1bSJed Brown     args: -dm_plex_hdf5_force_sequential
1783*c4762a1bSJed Brown     args: -distribute
1784*c4762a1bSJed Brown     test: # seq load, simple partitioner
1785*c4762a1bSJed Brown       suffix: 7_seq_hdf5_simple
1786*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1787*c4762a1bSJed Brown       args: -interpolate none
1788*c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1789*c4762a1bSJed Brown       suffix: 7_seq_hdf5_int_simple
1790*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1791*c4762a1bSJed Brown       args: -interpolate after_create
1792*c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1793*c4762a1bSJed Brown       nsize: {{2 4 5}}
1794*c4762a1bSJed Brown       suffix: 7_seq_hdf5_int_metis
1795*c4762a1bSJed Brown       requires: parmetis
1796*c4762a1bSJed Brown       args: -interpolate after_create
1797*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1798*c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1799*c4762a1bSJed Brown       suffix: 7_seq_hdf5_simple_int
1800*c4762a1bSJed Brown       nsize: {{2 4 5}}
1801*c4762a1bSJed Brown       args: -interpolate after_distribute
1802*c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1803*c4762a1bSJed Brown       nsize: {{2 4 5}}
1804*c4762a1bSJed Brown       suffix: 7_seq_hdf5_metis_int
1805*c4762a1bSJed Brown       requires: parmetis
1806*c4762a1bSJed Brown       args: -interpolate after_distribute
1807*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1808*c4762a1bSJed Brown 
1809*c4762a1bSJed Brown   testset: # 7 HDF5 PARALLEL LOAD
1810*c4762a1bSJed Brown     requires: hdf5 !complex
1811*c4762a1bSJed Brown     nsize: {{2 4 5}}
1812*c4762a1bSJed Brown     args: -dm_plex_check_all
1813*c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1814*c4762a1bSJed Brown     test: # par load
1815*c4762a1bSJed Brown       suffix: 7_par_hdf5
1816*c4762a1bSJed Brown       args: -interpolate none
1817*c4762a1bSJed Brown     test: # par load, par interpolation
1818*c4762a1bSJed Brown       suffix: 7_par_hdf5_int
1819*c4762a1bSJed Brown       args: -interpolate after_create
1820*c4762a1bSJed Brown     test: # par load, parmetis repartitioner
1821*c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1822*c4762a1bSJed Brown       suffix: 7_par_hdf5_parmetis
1823*c4762a1bSJed Brown       requires: parmetis
1824*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1825*c4762a1bSJed Brown       args: -interpolate none
1826*c4762a1bSJed Brown     test: # par load, par interpolation, parmetis repartitioner
1827*c4762a1bSJed Brown       suffix: 7_par_hdf5_int_parmetis
1828*c4762a1bSJed Brown       requires: parmetis
1829*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1830*c4762a1bSJed Brown       args: -interpolate after_create
1831*c4762a1bSJed Brown     test: # par load, parmetis partitioner, par interpolation
1832*c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1833*c4762a1bSJed Brown       suffix: 7_par_hdf5_parmetis_int
1834*c4762a1bSJed Brown       requires: parmetis
1835*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1836*c4762a1bSJed Brown       args: -interpolate after_distribute
1837*c4762a1bSJed Brown 
1838*c4762a1bSJed Brown     test:
1839*c4762a1bSJed Brown       suffix: 7_hdf5_hierarch
1840*c4762a1bSJed Brown       requires: hdf5 ptscotch !complex
1841*c4762a1bSJed Brown       nsize: {{2 3 4}separate output}
1842*c4762a1bSJed Brown       args: -distribute
1843*c4762a1bSJed Brown       args: -interpolate after_create
1844*c4762a1bSJed Brown       args: -petscpartitioner_type matpartitioning -petscpartitioner_view ::ascii_info
1845*c4762a1bSJed Brown       args: -mat_partitioning_type hierarch -mat_partitioning_hierarchical_nfineparts 2
1846*c4762a1bSJed Brown       args: -mat_partitioning_hierarchical_coarseparttype ptscotch -mat_partitioning_hierarchical_fineparttype ptscotch
1847*c4762a1bSJed Brown 
1848*c4762a1bSJed Brown   test:
1849*c4762a1bSJed Brown     suffix: 8
1850*c4762a1bSJed Brown     requires: hdf5 !complex
1851*c4762a1bSJed Brown     nsize: 4
1852*c4762a1bSJed Brown     args: -filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.h5 -dm_plex_create_from_hdf5_xdmf
1853*c4762a1bSJed Brown     args: -distribute 0 -interpolate after_create
1854*c4762a1bSJed 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
1855*c4762a1bSJed Brown     args: -dm_plex_check_all
1856*c4762a1bSJed Brown     args: -custom_view
1857*c4762a1bSJed Brown 
1858*c4762a1bSJed Brown   testset: # 9 HDF5 SEQUANTIAL LOAD
1859*c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1860*c4762a1bSJed Brown     args: -dm_plex_check_all
1861*c4762a1bSJed 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
1862*c4762a1bSJed Brown     args: -dm_plex_hdf5_force_sequential
1863*c4762a1bSJed Brown     args: -distribute
1864*c4762a1bSJed Brown     test: # seq load, simple partitioner
1865*c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple
1866*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1867*c4762a1bSJed Brown       args: -interpolate none
1868*c4762a1bSJed Brown     test: # seq load, seq interpolation, simple partitioner
1869*c4762a1bSJed Brown       suffix: 9_seq_hdf5_int_simple
1870*c4762a1bSJed Brown       nsize: {{1 2 4 5}}
1871*c4762a1bSJed Brown       args: -interpolate after_create
1872*c4762a1bSJed Brown     test: # seq load, seq interpolation, metis partitioner
1873*c4762a1bSJed Brown       nsize: {{2 4 5}}
1874*c4762a1bSJed Brown       suffix: 9_seq_hdf5_int_metis
1875*c4762a1bSJed Brown       requires: parmetis
1876*c4762a1bSJed Brown       args: -interpolate after_create
1877*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1878*c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1879*c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple_int
1880*c4762a1bSJed Brown       nsize: {{2 4 5}}
1881*c4762a1bSJed Brown       args: -interpolate after_distribute
1882*c4762a1bSJed Brown     test: # seq load, simple partitioner, par interpolation
1883*c4762a1bSJed Brown       # This is like 9_seq_hdf5_simple_int but testing error output of DMPlexCheckPointSFHeavy().
1884*c4762a1bSJed Brown       # Once 9_seq_hdf5_simple_int gets fixed, this one gets broken.
1885*c4762a1bSJed Brown       # We can then provide an intentionally broken mesh instead.
1886*c4762a1bSJed Brown       TODO: This test is broken because PointSF is fixed.
1887*c4762a1bSJed Brown       suffix: 9_seq_hdf5_simple_int_err
1888*c4762a1bSJed Brown       nsize: 4
1889*c4762a1bSJed Brown       args: -interpolate after_distribute
1890*c4762a1bSJed Brown       filter: sed -e "/PETSC ERROR/,$$d"
1891*c4762a1bSJed Brown     test: # seq load, metis partitioner, par interpolation
1892*c4762a1bSJed Brown       nsize: {{2 4 5}}
1893*c4762a1bSJed Brown       suffix: 9_seq_hdf5_metis_int
1894*c4762a1bSJed Brown       requires: parmetis
1895*c4762a1bSJed Brown       args: -interpolate after_distribute
1896*c4762a1bSJed Brown       args: -petscpartitioner_type parmetis
1897*c4762a1bSJed Brown 
1898*c4762a1bSJed Brown   testset: # 9 HDF5 PARALLEL LOAD
1899*c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1900*c4762a1bSJed Brown     nsize: {{2 4 5}}
1901*c4762a1bSJed Brown     args: -dm_plex_check_all
1902*c4762a1bSJed 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
1903*c4762a1bSJed Brown     test: # par load
1904*c4762a1bSJed Brown       suffix: 9_par_hdf5
1905*c4762a1bSJed Brown       args: -interpolate none
1906*c4762a1bSJed Brown     test: # par load, par interpolation
1907*c4762a1bSJed Brown       suffix: 9_par_hdf5_int
1908*c4762a1bSJed Brown       args: -interpolate after_create
1909*c4762a1bSJed Brown     test: # par load, parmetis repartitioner
1910*c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1911*c4762a1bSJed Brown       suffix: 9_par_hdf5_parmetis
1912*c4762a1bSJed Brown       requires: parmetis
1913*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1914*c4762a1bSJed Brown       args: -interpolate none
1915*c4762a1bSJed Brown     test: # par load, par interpolation, parmetis repartitioner
1916*c4762a1bSJed Brown       suffix: 9_par_hdf5_int_parmetis
1917*c4762a1bSJed Brown       requires: parmetis
1918*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1919*c4762a1bSJed Brown       args: -interpolate after_create
1920*c4762a1bSJed Brown     test: # par load, parmetis partitioner, par interpolation
1921*c4762a1bSJed Brown       TODO: Parallel partitioning of uninterpolated meshes not supported
1922*c4762a1bSJed Brown       suffix: 9_par_hdf5_parmetis_int
1923*c4762a1bSJed Brown       requires: parmetis
1924*c4762a1bSJed Brown       args: -distribute -petscpartitioner_type parmetis
1925*c4762a1bSJed Brown       args: -interpolate after_distribute
1926*c4762a1bSJed Brown 
1927*c4762a1bSJed Brown   testset: # 10 HDF5 PARALLEL LOAD
1928*c4762a1bSJed Brown     requires: hdf5 !complex datafilespath
1929*c4762a1bSJed Brown     nsize: {{2 4 7}}
1930*c4762a1bSJed Brown     args: -dm_plex_check_all
1931*c4762a1bSJed 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
1932*c4762a1bSJed Brown     test: # par load, par interpolation
1933*c4762a1bSJed Brown       suffix: 10_par_hdf5_int
1934*c4762a1bSJed Brown       args: -interpolate after_create
1935*c4762a1bSJed Brown TEST*/
1936