xref: /petsc/src/dm/impls/plex/tests/ex7.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Tests for mesh interpolation\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /* List of test meshes
6c4762a1bSJed Brown 
7c4762a1bSJed Brown Triangle
8c4762a1bSJed Brown --------
9c4762a1bSJed Brown Test 0:
10c4762a1bSJed Brown Two triangles sharing a face
11c4762a1bSJed Brown 
12c4762a1bSJed Brown         4
13c4762a1bSJed Brown       / | \
14c4762a1bSJed Brown      /  |  \
15c4762a1bSJed Brown     /   |   \
16c4762a1bSJed Brown    2  0 | 1  5
17c4762a1bSJed Brown     \   |   /
18c4762a1bSJed Brown      \  |  /
19c4762a1bSJed Brown       \ | /
20c4762a1bSJed Brown         3
21c4762a1bSJed Brown 
22c4762a1bSJed Brown should become
23c4762a1bSJed Brown 
24c4762a1bSJed Brown         4
25c4762a1bSJed Brown       / | \
26c4762a1bSJed Brown      8  |  9
27c4762a1bSJed Brown     /   |   \
28c4762a1bSJed Brown    2  0 7 1  5
29c4762a1bSJed Brown     \   |   /
30c4762a1bSJed Brown      6  |  10
31c4762a1bSJed Brown       \ | /
32c4762a1bSJed Brown         3
33c4762a1bSJed Brown 
34c4762a1bSJed Brown Tetrahedron
35c4762a1bSJed Brown -----------
36c4762a1bSJed Brown Test 0:
37c4762a1bSJed Brown Two tets sharing a face
38c4762a1bSJed Brown 
39c4762a1bSJed Brown  cell   5 _______    cell
40c4762a1bSJed Brown  0    / | \      \       1
41c4762a1bSJed Brown      /  |  \      \
42c4762a1bSJed Brown     /   |   \      \
43c4762a1bSJed Brown    2----|----4-----6
44c4762a1bSJed Brown     \   |   /      /
45c4762a1bSJed Brown      \  |  /     /
46c4762a1bSJed Brown       \ | /      /
47c4762a1bSJed Brown         3-------
48c4762a1bSJed Brown 
49c4762a1bSJed Brown should become
50c4762a1bSJed Brown 
51c4762a1bSJed Brown  cell   5 _______    cell
52c4762a1bSJed Brown  0    / | \      \       1
53c4762a1bSJed Brown      /  |  \      \
54c4762a1bSJed Brown    17   |   18 13  22
55c4762a1bSJed Brown    / 8 19 10 \      \
56c4762a1bSJed Brown   /     |     \      \
57c4762a1bSJed Brown  2---14-|------4--20--6
58c4762a1bSJed Brown   \     |     /      /
59c4762a1bSJed Brown    \ 9  | 7  /      /
60c4762a1bSJed Brown    16   |   15 11  21
61c4762a1bSJed Brown      \  |  /      /
62c4762a1bSJed Brown       \ | /      /
63c4762a1bSJed Brown         3-------
64c4762a1bSJed Brown 
65c4762a1bSJed Brown Quadrilateral
66c4762a1bSJed Brown -------------
67c4762a1bSJed Brown Test 0:
68c4762a1bSJed Brown Two quads sharing a face
69c4762a1bSJed Brown 
70c4762a1bSJed Brown    5-------4-------7
71c4762a1bSJed Brown    |       |       |
72c4762a1bSJed Brown    |   0   |   1   |
73c4762a1bSJed Brown    |       |       |
74c4762a1bSJed Brown    2-------3-------6
75c4762a1bSJed Brown 
76c4762a1bSJed Brown should become
77c4762a1bSJed Brown 
78c4762a1bSJed Brown    5--10---4--14---7
79c4762a1bSJed Brown    |       |       |
80c4762a1bSJed Brown   11   0   9   1  13
81c4762a1bSJed Brown    |       |       |
82c4762a1bSJed Brown    2---8---3--12---6
83c4762a1bSJed Brown 
84c4762a1bSJed Brown Test 1:
85c4762a1bSJed Brown A quad and a triangle sharing a face
86c4762a1bSJed Brown 
87c4762a1bSJed Brown    5-------4
88c4762a1bSJed Brown    |       | \
89c4762a1bSJed Brown    |   0   |  \
90c4762a1bSJed Brown    |       | 1 \
91c4762a1bSJed Brown    2-------3----6
92c4762a1bSJed Brown 
93c4762a1bSJed Brown should become
94c4762a1bSJed Brown 
95c4762a1bSJed Brown    5---9---4
96c4762a1bSJed Brown    |       | \
97c4762a1bSJed Brown   10   0   8  12
98c4762a1bSJed Brown    |       | 1 \
99c4762a1bSJed Brown    2---7---3-11-6
100c4762a1bSJed Brown 
101c4762a1bSJed Brown Hexahedron
102c4762a1bSJed Brown ----------
103c4762a1bSJed Brown Test 0:
104c4762a1bSJed Brown Two hexes sharing a face
105c4762a1bSJed Brown 
106c4762a1bSJed Brown cell   9-------------8-------------13 cell
107c4762a1bSJed Brown 0     /|            /|            /|     1
108c4762a1bSJed Brown      / |   15      / |   21      / |
109c4762a1bSJed Brown     /  |          /  |          /  |
110c4762a1bSJed Brown    6-------------7-------------12  |
111c4762a1bSJed Brown    |   |     18  |   |     24  |   |
112c4762a1bSJed Brown    |   |         |   |         |   |
113c4762a1bSJed Brown    |19 |         |17 |         |23 |
114c4762a1bSJed Brown    |   |  16     |   |   22    |   |
115c4762a1bSJed Brown    |   5---------|---4---------|---11
116c4762a1bSJed Brown    |  /          |  /          |  /
117c4762a1bSJed Brown    | /   14      | /    20     | /
118c4762a1bSJed Brown    |/            |/            |/
119c4762a1bSJed Brown    2-------------3-------------10
120c4762a1bSJed Brown 
121c4762a1bSJed Brown should become
122c4762a1bSJed Brown 
123c4762a1bSJed Brown cell   9-----31------8-----42------13 cell
124c4762a1bSJed Brown 0     /|            /|            /|     1
125c4762a1bSJed Brown     32 |   15      30|   21      41|
126c4762a1bSJed Brown     /  |          /  |          /  |
127c4762a1bSJed Brown    6-----29------7-----40------12  |
128c4762a1bSJed Brown    |   |     17  |   |     23  |   |
129c4762a1bSJed Brown    |  35         |  36         |   44
130c4762a1bSJed Brown    |19 |         |18 |         |24 |
131c4762a1bSJed Brown   34   |  16    33   |   22   43   |
132c4762a1bSJed Brown    |   5-----26--|---4-----37--|---11
133c4762a1bSJed Brown    |  /          |  /          |  /
134c4762a1bSJed Brown    | 25   14     | 27    20    | 38
135c4762a1bSJed Brown    |/            |/            |/
136c4762a1bSJed Brown    2-----28------3-----39------10
137c4762a1bSJed Brown 
138c4762a1bSJed Brown */
139c4762a1bSJed Brown 
140c4762a1bSJed Brown typedef struct {
141c4762a1bSJed Brown   PetscInt  testNum;      /* Indicates the mesh to create */
142c4762a1bSJed Brown   PetscInt  dim;          /* The topological mesh dimension */
14330602db0SMatthew G. Knepley   PetscBool simplex;      /* Use simplices or hexes */
144c4762a1bSJed Brown   PetscBool useGenerator; /* Construct mesh with a mesh generator */
145c4762a1bSJed Brown } AppCtx;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
148c4762a1bSJed Brown {
149c4762a1bSJed Brown   PetscFunctionBegin;
150c4762a1bSJed Brown   options->testNum      = 0;
151c4762a1bSJed Brown   options->dim          = 2;
15230602db0SMatthew G. Knepley   options->simplex      = PETSC_TRUE;
153c4762a1bSJed Brown   options->useGenerator = PETSC_FALSE;
154c4762a1bSJed Brown 
155d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
1569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL,0));
1579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL,1,3));
1589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL));
1599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL));
160d0609cedSBarry Smith   PetscOptionsEnd();
161c4762a1bSJed Brown   PetscFunctionReturn(0);
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
164c4762a1bSJed Brown PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM dm)
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
167c4762a1bSJed Brown   PetscMPIInt    rank;
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   PetscFunctionBegin;
1709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
171dd400576SPatrick Sanan   if (rank == 0) {
172c4762a1bSJed Brown     switch (testNum) {
173c4762a1bSJed Brown     case 0:
174c4762a1bSJed Brown     {
175c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
176c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
177c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
178c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
179c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
180c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
183c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
1849566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
185c4762a1bSJed Brown       }
186c4762a1bSJed Brown     }
187c4762a1bSJed Brown     break;
188c4762a1bSJed Brown     default:
18963a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown   } else {
192c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
1959566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown   PetscFunctionReturn(0);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
203c4762a1bSJed Brown   PetscMPIInt    rank;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   PetscFunctionBegin;
2069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
207dd400576SPatrick Sanan   if (rank == 0) {
208c4762a1bSJed Brown     switch (testNum) {
209c4762a1bSJed Brown     case 0:
210c4762a1bSJed Brown     {
211c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
212c4762a1bSJed Brown       PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
213c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 4, 3, 5,  3, 4, 6, 5};
214c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
215c4762a1bSJed Brown       PetscScalar vertexCoords[15]    = {0.0, 0.0, -0.5,  0.0, -0.5, 0.0,  1.0, 0.0, 0.0,  0.0, 0.5, 0.0,  0.0, 0.0, 0.5};
216c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
217c4762a1bSJed Brown 
2189566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
219c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
2209566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
221c4762a1bSJed Brown       }
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown     break;
224c4762a1bSJed Brown     default:
22563a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
226c4762a1bSJed Brown     }
227c4762a1bSJed Brown   } else {
228c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
229c4762a1bSJed Brown 
2309566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2319566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
232c4762a1bSJed Brown   }
233c4762a1bSJed Brown   PetscFunctionReturn(0);
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
236c4762a1bSJed Brown PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm)
237c4762a1bSJed Brown {
238c4762a1bSJed Brown   PetscInt       depth = 1, p;
239c4762a1bSJed Brown   PetscMPIInt    rank;
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   PetscFunctionBegin;
2429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
243dd400576SPatrick Sanan   if (rank == 0) {
244c4762a1bSJed Brown     switch (testNum) {
245c4762a1bSJed Brown     case 0:
246c4762a1bSJed Brown     {
247c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
248c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
249c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
250c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
251c4762a1bSJed Brown       PetscScalar vertexCoords[12]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0, 0.5, 1.0};
252c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
253c4762a1bSJed Brown 
2549566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
255c4762a1bSJed Brown       for (p = 0; p < 6; ++p) {
2569566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
257c4762a1bSJed Brown       }
258c4762a1bSJed Brown     }
259c4762a1bSJed Brown     break;
260c4762a1bSJed Brown     case 1:
261c4762a1bSJed Brown     {
262c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
263c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 3, 0, 0, 0, 0, 0};
264c4762a1bSJed Brown       PetscInt    cones[7]            = {2, 3, 4, 5,  3, 6, 4};
265c4762a1bSJed Brown       PetscInt    coneOrientations[7] = {0, 0, 0, 0,  0, 0, 0};
266c4762a1bSJed Brown       PetscScalar vertexCoords[10]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0};
267c4762a1bSJed Brown       PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
268c4762a1bSJed Brown 
2699566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
270c4762a1bSJed Brown       for (p = 0; p < 5; ++p) {
2719566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
272c4762a1bSJed Brown       }
273c4762a1bSJed Brown     }
274c4762a1bSJed Brown     break;
275c4762a1bSJed Brown     default:
27663a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
277c4762a1bSJed Brown     }
278c4762a1bSJed Brown   } else {
279c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
280c4762a1bSJed Brown 
2819566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2829566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown   PetscFunctionReturn(0);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm)
288c4762a1bSJed Brown {
289c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
290c4762a1bSJed Brown   PetscMPIInt    rank;
291c4762a1bSJed Brown 
292c4762a1bSJed Brown   PetscFunctionBegin;
2939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
294dd400576SPatrick Sanan   if (rank == 0) {
295c4762a1bSJed Brown     switch (testNum) {
296c4762a1bSJed Brown     case 0:
297c4762a1bSJed Brown     {
298c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
299c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0};
300c4762a1bSJed Brown       PetscInt    cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
301c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0, 0,0, 0, 0,0};
302c4762a1bSJed Brown       PetscScalar vertexCoords[36]     = {-0.5,0.0,0.0, 0.0,0.0,0.0, 0.0,1.0,0.0, -0.5,1.0,0.0,
303c4762a1bSJed Brown                                           -0.5,0.0,1.0, 0.0,0.0,1.0, 0.0,1.0,1.0, -0.5,1.0,1.0,
304c4762a1bSJed Brown                                            0.5,0.0,0.0, 0.5,1.0,0.0, 0.5,0.0,1.0,  0.5,1.0,1.0};
305c4762a1bSJed Brown       PetscInt    markerPoints[16]     = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
306c4762a1bSJed Brown 
3079566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
308c4762a1bSJed Brown       for (p = 0; p < 8; ++p) {
3099566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
310c4762a1bSJed Brown       }
311c4762a1bSJed Brown     }
312c4762a1bSJed Brown     break;
313c4762a1bSJed Brown     default:
31463a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
315c4762a1bSJed Brown     }
316c4762a1bSJed Brown   } else {
317c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
318c4762a1bSJed Brown 
3199566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
3209566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
321c4762a1bSJed Brown   }
322c4762a1bSJed Brown   PetscFunctionReturn(0);
323c4762a1bSJed Brown }
324c4762a1bSJed Brown 
325c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm)
326c4762a1bSJed Brown {
327c4762a1bSJed Brown   PetscBool      useGenerator = user->useGenerator;
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
3319566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
332b5a892a1SMatthew G. Knepley   if (!useGenerator) {
33330602db0SMatthew G. Knepley     PetscInt  dim     = user->dim;
33430602db0SMatthew G. Knepley     PetscBool simplex = user->simplex;
33530602db0SMatthew G. Knepley 
3369566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dm, dim));
337c4762a1bSJed Brown     switch (dim) {
338c4762a1bSJed Brown     case 2:
33930602db0SMatthew G. Knepley       if (simplex) {
3409566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_2D(comm, *dm));
341c4762a1bSJed Brown       } else {
3429566063dSJacob Faibussowitsch         PetscCall(CreateQuad_2D(comm, testNum, *dm));
343c4762a1bSJed Brown       }
344c4762a1bSJed Brown       break;
345c4762a1bSJed Brown     case 3:
34630602db0SMatthew G. Knepley       if (simplex) {
3479566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_3D(comm, *dm));
348c4762a1bSJed Brown       } else {
3499566063dSJacob Faibussowitsch         PetscCall(CreateHex_3D(comm, *dm));
350c4762a1bSJed Brown       }
351c4762a1bSJed Brown       break;
352c4762a1bSJed Brown     default:
35363a3b9bcSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
354c4762a1bSJed Brown     }
355c4762a1bSJed Brown   }
3569566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
3579566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, "Interpolated Mesh"));
3589566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
359c4762a1bSJed Brown   PetscFunctionReturn(0);
360c4762a1bSJed Brown }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown int main(int argc, char **argv)
363c4762a1bSJed Brown {
36430602db0SMatthew G. Knepley   DM             dm;
365c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
366c4762a1bSJed Brown 
367*327415f7SBarry Smith   PetscFunctionBeginUser;
3689566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
3699566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3709566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm));
3719566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3729566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
373b122ec5aSJacob Faibussowitsch   return 0;
374c4762a1bSJed Brown }
375c4762a1bSJed Brown 
376c4762a1bSJed Brown /*TEST
377b5a892a1SMatthew G. Knepley   testset:
378b5a892a1SMatthew G. Knepley     args: -dm_plex_interpolate_pre -dm_plex_check_all
379c4762a1bSJed Brown 
380c4762a1bSJed Brown     # Two cell test meshes 0-7
381c4762a1bSJed Brown     test:
382c4762a1bSJed Brown       suffix: 0
38367d101c4SVaclav Hapla       args: -dim 2 -dm_view ascii::ascii_info_detail
384c4762a1bSJed Brown     test:
385c4762a1bSJed Brown       suffix: 1
386c4762a1bSJed Brown       nsize: 2
387e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail
388c4762a1bSJed Brown     test:
389c4762a1bSJed Brown       suffix: 2
39030602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
391c4762a1bSJed Brown     test:
392c4762a1bSJed Brown       suffix: 3
393c4762a1bSJed Brown       nsize: 2
394e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
395c4762a1bSJed Brown     test:
396c4762a1bSJed Brown       suffix: 4
397c4762a1bSJed Brown       args: -dim 3 -dm_view ascii::ascii_info_detail
398c4762a1bSJed Brown     test:
399c4762a1bSJed Brown       suffix: 5
400c4762a1bSJed Brown       nsize: 2
401e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail
402c4762a1bSJed Brown     test:
403c4762a1bSJed Brown       suffix: 6
40430602db0SMatthew G. Knepley       args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
405c4762a1bSJed Brown     test:
406c4762a1bSJed Brown       suffix: 7
407c4762a1bSJed Brown       nsize: 2
408e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
409c4762a1bSJed Brown     # 2D Hybrid Mesh 8
410c4762a1bSJed Brown     test:
411c4762a1bSJed Brown       suffix: 8
41230602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail
413b5a892a1SMatthew G. Knepley 
414b5a892a1SMatthew G. Knepley   testset:
415b5a892a1SMatthew G. Knepley     args: -dm_plex_check_all
416b5a892a1SMatthew G. Knepley 
417da9060c4SMatthew G. Knepley     # Reference cells
418da9060c4SMatthew G. Knepley     test:
419da9060c4SMatthew G. Knepley       suffix: 12
420b5a892a1SMatthew G. Knepley       args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all
421c4762a1bSJed Brown     # TetGen meshes 9-10
422c4762a1bSJed Brown     test:
423c4762a1bSJed Brown       suffix: 9
424c4762a1bSJed Brown       requires: triangle
42530602db0SMatthew G. Knepley       args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0
426c4762a1bSJed Brown     test:
427c4762a1bSJed Brown       suffix: 10
428c4762a1bSJed Brown       requires: ctetgen
42930602db0SMatthew G. Knepley       args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0
430c4762a1bSJed Brown     # Cubit meshes 11
431c4762a1bSJed Brown     test:
432c4762a1bSJed Brown       suffix: 11
433c4762a1bSJed Brown       requires: exodusii
43430602db0SMatthew G. Knepley       args: -use_generator -dm_plex_filename ${wPETSC_DIR}/share/petsc/datafiles/meshes/blockcylinder-50.exo -dm_view ascii::ascii_info_detail -dm_coord_space 0
435c4762a1bSJed Brown 
436c4762a1bSJed Brown TEST*/
437