xref: /petsc/src/dm/impls/plex/tests/ex7.c (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
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);
158c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
16030602db0SMatthew G. Knepley   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
161c4762a1bSJed Brown   ierr = PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL);CHKERRQ(ierr);
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   PetscErrorCode ierr;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBegin;
173ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
174dd400576SPatrick Sanan   if (rank == 0) {
175c4762a1bSJed Brown     switch (testNum) {
176c4762a1bSJed Brown     case 0:
177c4762a1bSJed Brown     {
178c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
179c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
180c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
181c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
182c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
183c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
184c4762a1bSJed Brown 
185c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
186c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
187c4762a1bSJed Brown         ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
188c4762a1bSJed Brown       }
189c4762a1bSJed Brown     }
190c4762a1bSJed Brown     break;
191c4762a1bSJed Brown     default:
19298921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
193c4762a1bSJed Brown     }
194c4762a1bSJed Brown   } else {
195c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
196c4762a1bSJed Brown 
197c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
198b5a892a1SMatthew G. Knepley     ierr = DMCreateLabel(dm, "marker");CHKERRQ(ierr);
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown   PetscFunctionReturn(0);
201c4762a1bSJed Brown }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm)
204c4762a1bSJed Brown {
205c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
206c4762a1bSJed Brown   PetscMPIInt    rank;
207c4762a1bSJed Brown   PetscErrorCode ierr;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
210ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
211dd400576SPatrick Sanan   if (rank == 0) {
212c4762a1bSJed Brown     switch (testNum) {
213c4762a1bSJed Brown     case 0:
214c4762a1bSJed Brown     {
215c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
216c4762a1bSJed Brown       PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
217c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 4, 3, 5,  3, 4, 6, 5};
218c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
219c4762a1bSJed 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};
220c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
221c4762a1bSJed Brown 
222c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
223c4762a1bSJed Brown       for (p = 0; p < 4; ++p) {
224c4762a1bSJed Brown         ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
225c4762a1bSJed Brown       }
226c4762a1bSJed Brown     }
227c4762a1bSJed Brown     break;
228c4762a1bSJed Brown     default:
22998921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
230c4762a1bSJed Brown     }
231c4762a1bSJed Brown   } else {
232c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
233c4762a1bSJed Brown 
234c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
235b5a892a1SMatthew G. Knepley     ierr = DMCreateLabel(dm, "marker");CHKERRQ(ierr);
236c4762a1bSJed Brown   }
237c4762a1bSJed Brown   PetscFunctionReturn(0);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   PetscInt       depth = 1, p;
243c4762a1bSJed Brown   PetscMPIInt    rank;
244c4762a1bSJed Brown   PetscErrorCode ierr;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   PetscFunctionBegin;
247ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
248dd400576SPatrick Sanan   if (rank == 0) {
249c4762a1bSJed Brown     switch (testNum) {
250c4762a1bSJed Brown     case 0:
251c4762a1bSJed Brown     {
252c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
253c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
254c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4};
255c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
256c4762a1bSJed 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};
257c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
258c4762a1bSJed Brown 
259c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
260c4762a1bSJed Brown       for (p = 0; p < 6; ++p) {
261c4762a1bSJed Brown         ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
262c4762a1bSJed Brown       }
263c4762a1bSJed Brown     }
264c4762a1bSJed Brown     break;
265c4762a1bSJed Brown     case 1:
266c4762a1bSJed Brown     {
267c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
268c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 3, 0, 0, 0, 0, 0};
269c4762a1bSJed Brown       PetscInt    cones[7]            = {2, 3, 4, 5,  3, 6, 4};
270c4762a1bSJed Brown       PetscInt    coneOrientations[7] = {0, 0, 0, 0,  0, 0, 0};
271c4762a1bSJed Brown       PetscScalar vertexCoords[10]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0};
272c4762a1bSJed Brown       PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
273c4762a1bSJed Brown 
274c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
275c4762a1bSJed Brown       for (p = 0; p < 5; ++p) {
276c4762a1bSJed Brown         ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
277c4762a1bSJed Brown       }
278c4762a1bSJed Brown     }
279c4762a1bSJed Brown     break;
280c4762a1bSJed Brown     default:
28198921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown   } else {
284c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
285c4762a1bSJed Brown 
286c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
287b5a892a1SMatthew G. Knepley     ierr = DMCreateLabel(dm, "marker");CHKERRQ(ierr);
288c4762a1bSJed Brown   }
289c4762a1bSJed Brown   PetscFunctionReturn(0);
290c4762a1bSJed Brown }
291c4762a1bSJed Brown 
292c4762a1bSJed Brown PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm)
293c4762a1bSJed Brown {
294c4762a1bSJed Brown   PetscInt       depth = 1, testNum  = 0, p;
295c4762a1bSJed Brown   PetscMPIInt    rank;
296c4762a1bSJed Brown   PetscErrorCode ierr;
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   PetscFunctionBegin;
299ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
300dd400576SPatrick Sanan   if (rank == 0) {
301c4762a1bSJed Brown     switch (testNum) {
302c4762a1bSJed Brown     case 0:
303c4762a1bSJed Brown     {
304c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
305c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0,0,0,0,0,0,0,0,0,0,0,0};
306c4762a1bSJed Brown       PetscInt    cones[16]            = {2,5,4,3,6,7,8,9,  3,4,11,10,7,12,13,8};
307c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0,0,0,0,0,0,0,0,  0,0, 0, 0,0, 0, 0,0};
308c4762a1bSJed 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,
309c4762a1bSJed 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,
310c4762a1bSJed 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};
311c4762a1bSJed Brown       PetscInt    markerPoints[16]     = {2,1,3,1,4,1,5,1,6,1,7,1,8,1,9,1};
312c4762a1bSJed Brown 
313c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
314c4762a1bSJed Brown       for (p = 0; p < 8; ++p) {
315c4762a1bSJed Brown         ierr = DMSetLabelValue(dm, "marker", markerPoints[p*2], markerPoints[p*2+1]);CHKERRQ(ierr);
316c4762a1bSJed Brown       }
317c4762a1bSJed Brown     }
318c4762a1bSJed Brown     break;
319c4762a1bSJed Brown     default:
32098921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   } else {
323c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
326b5a892a1SMatthew G. Knepley     ierr = DMCreateLabel(dm, "marker");CHKERRQ(ierr);
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown   PetscFunctionReturn(0);
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm)
332c4762a1bSJed Brown {
333c4762a1bSJed Brown   PetscBool      useGenerator = user->useGenerator;
334c4762a1bSJed Brown   PetscErrorCode ierr;
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   PetscFunctionBegin;
337c4762a1bSJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
338c4762a1bSJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
339b5a892a1SMatthew G. Knepley   if (!useGenerator) {
34030602db0SMatthew G. Knepley     PetscInt  dim     = user->dim;
34130602db0SMatthew G. Knepley     PetscBool simplex = user->simplex;
34230602db0SMatthew G. Knepley 
343c4762a1bSJed Brown     ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
344c4762a1bSJed Brown     switch (dim) {
345c4762a1bSJed Brown     case 2:
34630602db0SMatthew G. Knepley       if (simplex) {
347c4762a1bSJed Brown         ierr = CreateSimplex_2D(comm, *dm);CHKERRQ(ierr);
348c4762a1bSJed Brown       } else {
349c4762a1bSJed Brown         ierr = CreateQuad_2D(comm, testNum, *dm);CHKERRQ(ierr);
350c4762a1bSJed Brown       }
351c4762a1bSJed Brown       break;
352c4762a1bSJed Brown     case 3:
35330602db0SMatthew G. Knepley       if (simplex) {
354c4762a1bSJed Brown         ierr = CreateSimplex_3D(comm, *dm);CHKERRQ(ierr);
355c4762a1bSJed Brown       } else {
356c4762a1bSJed Brown         ierr = CreateHex_3D(comm, *dm);CHKERRQ(ierr);
357c4762a1bSJed Brown       }
358c4762a1bSJed Brown       break;
359c4762a1bSJed Brown     default:
36098921bdaSJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
361c4762a1bSJed Brown     }
362c4762a1bSJed Brown   }
363da9060c4SMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
364b5a892a1SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dm, "Interpolated Mesh");CHKERRQ(ierr);
365c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
366c4762a1bSJed Brown   PetscFunctionReturn(0);
367c4762a1bSJed Brown }
368c4762a1bSJed Brown 
369c4762a1bSJed Brown int main(int argc, char **argv)
370c4762a1bSJed Brown {
37130602db0SMatthew G. Knepley   DM             dm;
372c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
373c4762a1bSJed Brown   PetscErrorCode ierr;
374c4762a1bSJed Brown 
375c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
376c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
37730602db0SMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm);CHKERRQ(ierr);
37830602db0SMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
379c4762a1bSJed Brown   ierr = PetscFinalize();
380c4762a1bSJed Brown   return ierr;
381c4762a1bSJed Brown }
382c4762a1bSJed Brown 
383c4762a1bSJed Brown /*TEST
384b5a892a1SMatthew G. Knepley   testset:
385b5a892a1SMatthew G. Knepley     args: -dm_plex_interpolate_pre -dm_plex_check_all
386c4762a1bSJed Brown 
387c4762a1bSJed Brown     # Two cell test meshes 0-7
388c4762a1bSJed Brown     test:
389c4762a1bSJed Brown       suffix: 0
39067d101c4SVaclav Hapla       args: -dim 2 -dm_view ascii::ascii_info_detail
391c4762a1bSJed Brown     test:
392c4762a1bSJed Brown       suffix: 1
393c4762a1bSJed Brown       nsize: 2
394*e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail
395c4762a1bSJed Brown     test:
396c4762a1bSJed Brown       suffix: 2
39730602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
398c4762a1bSJed Brown     test:
399c4762a1bSJed Brown       suffix: 3
400c4762a1bSJed Brown       nsize: 2
401*e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
402c4762a1bSJed Brown     test:
403c4762a1bSJed Brown       suffix: 4
404c4762a1bSJed Brown       args: -dim 3 -dm_view ascii::ascii_info_detail
405c4762a1bSJed Brown     test:
406c4762a1bSJed Brown       suffix: 5
407c4762a1bSJed Brown       nsize: 2
408*e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail
409c4762a1bSJed Brown     test:
410c4762a1bSJed Brown       suffix: 6
41130602db0SMatthew G. Knepley       args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
412c4762a1bSJed Brown     test:
413c4762a1bSJed Brown       suffix: 7
414c4762a1bSJed Brown       nsize: 2
415*e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
416c4762a1bSJed Brown     # 2D Hybrid Mesh 8
417c4762a1bSJed Brown     test:
418c4762a1bSJed Brown       suffix: 8
41930602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail
420b5a892a1SMatthew G. Knepley 
421b5a892a1SMatthew G. Knepley   testset:
422b5a892a1SMatthew G. Knepley     args: -dm_plex_check_all
423b5a892a1SMatthew G. Knepley 
424da9060c4SMatthew G. Knepley     # Reference cells
425da9060c4SMatthew G. Knepley     test:
426da9060c4SMatthew G. Knepley       suffix: 12
427b5a892a1SMatthew G. Knepley       args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all
428c4762a1bSJed Brown     # TetGen meshes 9-10
429c4762a1bSJed Brown     test:
430c4762a1bSJed Brown       suffix: 9
431c4762a1bSJed Brown       requires: triangle
43230602db0SMatthew G. Knepley       args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0
433c4762a1bSJed Brown     test:
434c4762a1bSJed Brown       suffix: 10
435c4762a1bSJed Brown       requires: ctetgen
43630602db0SMatthew G. Knepley       args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0
437c4762a1bSJed Brown     # Cubit meshes 11
438c4762a1bSJed Brown     test:
439c4762a1bSJed Brown       suffix: 11
440c4762a1bSJed Brown       requires: exodusii
44130602db0SMatthew 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
442c4762a1bSJed Brown 
443c4762a1bSJed Brown TEST*/
444