xref: /petsc/src/dm/impls/plex/tests/ex7.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscErrorCode ierr;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   PetscFunctionBegin;
152c4762a1bSJed Brown   options->testNum      = 0;
153c4762a1bSJed Brown   options->dim          = 2;
15430602db0SMatthew G. Knepley   options->simplex      = PETSC_TRUE;
155c4762a1bSJed Brown   options->useGenerator = PETSC_FALSE;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");CHKERRQ(ierr);
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL,0));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL,1,3));
1605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL));
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL));
1621e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM dm)
167c4762a1bSJed Brown {
168c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
169c4762a1bSJed Brown   PetscMPIInt    rank;
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   PetscFunctionBegin;
1725f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
173dd400576SPatrick Sanan   if (rank == 0) {
174c4762a1bSJed Brown     switch (testNum) {
175c4762a1bSJed Brown     case 0:
176c4762a1bSJed Brown     {
177c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
178c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
179c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
180c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
181c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
182c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
183c4762a1bSJed Brown 
1845f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
185c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
1865f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
187c4762a1bSJed Brown       }
188c4762a1bSJed Brown     }
189c4762a1bSJed Brown     break;
190c4762a1bSJed Brown     default:
19198921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
192c4762a1bSJed Brown     }
193c4762a1bSJed Brown   } else {
194c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
195c4762a1bSJed Brown 
1965f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
1975f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "marker"));
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown   PetscFunctionReturn(0);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm)
203c4762a1bSJed Brown {
204c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
205c4762a1bSJed Brown   PetscMPIInt    rank;
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   PetscFunctionBegin;
2085f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
209dd400576SPatrick Sanan   if (rank == 0) {
210c4762a1bSJed Brown     switch (testNum) {
211c4762a1bSJed Brown     case 0:
212c4762a1bSJed Brown     {
213c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
214c4762a1bSJed Brown       PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
215c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 4, 3, 5,  3, 4, 6, 5};
216c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
217c4762a1bSJed 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};
218c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
219c4762a1bSJed Brown 
2205f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
221c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
2225f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
223c4762a1bSJed Brown       }
224c4762a1bSJed Brown     }
225c4762a1bSJed Brown     break;
226c4762a1bSJed Brown     default:
22798921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
228c4762a1bSJed Brown     }
229c4762a1bSJed Brown   } else {
230c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
231c4762a1bSJed Brown 
2325f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2335f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "marker"));
234c4762a1bSJed Brown   }
235c4762a1bSJed Brown   PetscFunctionReturn(0);
236c4762a1bSJed Brown }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm)
239c4762a1bSJed Brown {
240c4762a1bSJed Brown   PetscInt       depth = 1, p;
241c4762a1bSJed Brown   PetscMPIInt    rank;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   PetscFunctionBegin;
2445f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
245dd400576SPatrick Sanan   if (rank == 0) {
246c4762a1bSJed Brown     switch (testNum) {
247c4762a1bSJed Brown     case 0:
248c4762a1bSJed Brown     {
249c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
250c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
251c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
252c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
253c4762a1bSJed 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};
254c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
255c4762a1bSJed Brown 
2565f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
257c4762a1bSJed Brown       for (p = 0; p < 6; ++p) {
2585f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
259c4762a1bSJed Brown       }
260c4762a1bSJed Brown     }
261c4762a1bSJed Brown     break;
262c4762a1bSJed Brown     case 1:
263c4762a1bSJed Brown     {
264c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
265c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 3, 0, 0, 0, 0, 0};
266c4762a1bSJed Brown       PetscInt    cones[7]            = {2, 3, 4, 5,  3, 6, 4};
267c4762a1bSJed Brown       PetscInt    coneOrientations[7] = {0, 0, 0, 0,  0, 0, 0};
268c4762a1bSJed Brown       PetscScalar vertexCoords[10]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0};
269c4762a1bSJed Brown       PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
270c4762a1bSJed Brown 
2715f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
272c4762a1bSJed Brown       for (p = 0; p < 5; ++p) {
2735f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
274c4762a1bSJed Brown       }
275c4762a1bSJed Brown     }
276c4762a1bSJed Brown     break;
277c4762a1bSJed Brown     default:
27898921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
279c4762a1bSJed Brown     }
280c4762a1bSJed Brown   } else {
281c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
282c4762a1bSJed Brown 
2835f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2845f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "marker"));
285c4762a1bSJed Brown   }
286c4762a1bSJed Brown   PetscFunctionReturn(0);
287c4762a1bSJed Brown }
288c4762a1bSJed Brown 
289c4762a1bSJed Brown PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm)
290c4762a1bSJed Brown {
291c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
292c4762a1bSJed Brown   PetscMPIInt    rank;
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   PetscFunctionBegin;
2955f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(comm, &rank));
296dd400576SPatrick Sanan   if (rank == 0) {
297c4762a1bSJed Brown     switch (testNum) {
298c4762a1bSJed Brown     case 0:
299c4762a1bSJed Brown     {
300c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
301c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0};
302c4762a1bSJed Brown       PetscInt    cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
303c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0, 0,0, 0, 0,0};
304c4762a1bSJed 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,
305c4762a1bSJed 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,
306c4762a1bSJed 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};
307c4762a1bSJed Brown       PetscInt    markerPoints[16]     = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
308c4762a1bSJed Brown 
3095f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
310c4762a1bSJed Brown       for (p = 0; p < 8; ++p) {
3115f80ce2aSJacob Faibussowitsch         CHKERRQ(DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]));
312c4762a1bSJed Brown       }
313c4762a1bSJed Brown     }
314c4762a1bSJed Brown     break;
315c4762a1bSJed Brown     default:
31698921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
317c4762a1bSJed Brown     }
318c4762a1bSJed Brown   } else {
319c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
320c4762a1bSJed Brown 
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
3225f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateLabel(dm, "marker"));
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown   PetscFunctionReturn(0);
325c4762a1bSJed Brown }
326c4762a1bSJed Brown 
327c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm)
328c4762a1bSJed Brown {
329c4762a1bSJed Brown   PetscBool      useGenerator = user->useGenerator;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   PetscFunctionBegin;
3325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm, dm));
3335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm, DMPLEX));
334b5a892a1SMatthew G. Knepley   if (!useGenerator) {
33530602db0SMatthew G. Knepley     PetscInt  dim     = user->dim;
33630602db0SMatthew G. Knepley     PetscBool simplex = user->simplex;
33730602db0SMatthew G. Knepley 
3385f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetDimension(*dm, dim));
339c4762a1bSJed Brown     switch (dim) {
340c4762a1bSJed Brown     case 2:
34130602db0SMatthew G. Knepley       if (simplex) {
3425f80ce2aSJacob Faibussowitsch         CHKERRQ(CreateSimplex_2D(comm, *dm));
343c4762a1bSJed Brown       } else {
3445f80ce2aSJacob Faibussowitsch         CHKERRQ(CreateQuad_2D(comm, testNum, *dm));
345c4762a1bSJed Brown       }
346c4762a1bSJed Brown       break;
347c4762a1bSJed Brown     case 3:
34830602db0SMatthew G. Knepley       if (simplex) {
3495f80ce2aSJacob Faibussowitsch         CHKERRQ(CreateSimplex_3D(comm, *dm));
350c4762a1bSJed Brown       } else {
3515f80ce2aSJacob Faibussowitsch         CHKERRQ(CreateHex_3D(comm, *dm));
352c4762a1bSJed Brown       }
353c4762a1bSJed Brown       break;
354c4762a1bSJed Brown     default:
35598921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
356c4762a1bSJed Brown     }
357c4762a1bSJed Brown   }
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(*dm));
3595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectSetName((PetscObject) *dm, "Interpolated Mesh"));
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(*dm, NULL, "-dm_view"));
361c4762a1bSJed Brown   PetscFunctionReturn(0);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
364c4762a1bSJed Brown int main(int argc, char **argv)
365c4762a1bSJed Brown {
36630602db0SMatthew G. Knepley   DM             dm;
367c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
368c4762a1bSJed Brown 
369*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
3705f80ce2aSJacob Faibussowitsch   CHKERRQ(ProcessOptions(PETSC_COMM_WORLD, &user));
3715f80ce2aSJacob Faibussowitsch   CHKERRQ(CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm));
3725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
373*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
374*b122ec5aSJacob Faibussowitsch   return 0;
375c4762a1bSJed Brown }
376c4762a1bSJed Brown 
377c4762a1bSJed Brown /*TEST
378b5a892a1SMatthew G. Knepley   testset:
379b5a892a1SMatthew G. Knepley     args: -dm_plex_interpolate_pre -dm_plex_check_all
380c4762a1bSJed Brown 
381c4762a1bSJed Brown     # Two cell test meshes 0-7
382c4762a1bSJed Brown     test:
383c4762a1bSJed Brown       suffix: 0
38467d101c4SVaclav Hapla       args: -dim 2 -dm_view ascii::ascii_info_detail
385c4762a1bSJed Brown     test:
386c4762a1bSJed Brown       suffix: 1
387c4762a1bSJed Brown       nsize: 2
388e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail
389c4762a1bSJed Brown     test:
390c4762a1bSJed Brown       suffix: 2
39130602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
392c4762a1bSJed Brown     test:
393c4762a1bSJed Brown       suffix: 3
394c4762a1bSJed Brown       nsize: 2
395e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
396c4762a1bSJed Brown     test:
397c4762a1bSJed Brown       suffix: 4
398c4762a1bSJed Brown       args: -dim 3 -dm_view ascii::ascii_info_detail
399c4762a1bSJed Brown     test:
400c4762a1bSJed Brown       suffix: 5
401c4762a1bSJed Brown       nsize: 2
402e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail
403c4762a1bSJed Brown     test:
404c4762a1bSJed Brown       suffix: 6
40530602db0SMatthew G. Knepley       args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
406c4762a1bSJed Brown     test:
407c4762a1bSJed Brown       suffix: 7
408c4762a1bSJed Brown       nsize: 2
409e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
410c4762a1bSJed Brown     # 2D Hybrid Mesh 8
411c4762a1bSJed Brown     test:
412c4762a1bSJed Brown       suffix: 8
41330602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail
414b5a892a1SMatthew G. Knepley 
415b5a892a1SMatthew G. Knepley   testset:
416b5a892a1SMatthew G. Knepley     args: -dm_plex_check_all
417b5a892a1SMatthew G. Knepley 
418da9060c4SMatthew G. Knepley     # Reference cells
419da9060c4SMatthew G. Knepley     test:
420da9060c4SMatthew G. Knepley       suffix: 12
421b5a892a1SMatthew G. Knepley       args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all
422c4762a1bSJed Brown     # TetGen meshes 9-10
423c4762a1bSJed Brown     test:
424c4762a1bSJed Brown       suffix: 9
425c4762a1bSJed Brown       requires: triangle
42630602db0SMatthew G. Knepley       args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0
427c4762a1bSJed Brown     test:
428c4762a1bSJed Brown       suffix: 10
429c4762a1bSJed Brown       requires: ctetgen
43030602db0SMatthew G. Knepley       args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0
431c4762a1bSJed Brown     # Cubit meshes 11
432c4762a1bSJed Brown     test:
433c4762a1bSJed Brown       suffix: 11
434c4762a1bSJed Brown       requires: exodusii
43530602db0SMatthew 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
436c4762a1bSJed Brown 
437c4762a1bSJed Brown TEST*/
438