xref: /petsc/src/dm/impls/plex/tests/ex7.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
1479371c9d4SSatish Balay PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
148c4762a1bSJed Brown   PetscFunctionBegin;
149c4762a1bSJed Brown   options->testNum      = 0;
150c4762a1bSJed Brown   options->dim          = 2;
15130602db0SMatthew G. Knepley   options->simplex      = PETSC_TRUE;
152c4762a1bSJed Brown   options->useGenerator = PETSC_FALSE;
153c4762a1bSJed Brown 
154d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
1559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL, 0));
1569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL, 1, 3));
1579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL));
1589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL));
159d0609cedSBarry Smith   PetscOptionsEnd();
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
1639371c9d4SSatish Balay PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM dm) {
164c4762a1bSJed Brown   PetscInt    depth = 1, testNum = 0, p;
165c4762a1bSJed Brown   PetscMPIInt rank;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
169dd400576SPatrick Sanan   if (rank == 0) {
170c4762a1bSJed Brown     switch (testNum) {
1719371c9d4SSatish Balay     case 0: {
172c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
173c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
174c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
175c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
176c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
177c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
178c4762a1bSJed Brown 
1799566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
180*48a46eb9SPierre Jolivet       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
1819371c9d4SSatish Balay     } break;
1829371c9d4SSatish Balay     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown   } else {
185c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
1889566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
189c4762a1bSJed Brown   }
190c4762a1bSJed Brown   PetscFunctionReturn(0);
191c4762a1bSJed Brown }
192c4762a1bSJed Brown 
1939371c9d4SSatish Balay PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm) {
194c4762a1bSJed Brown   PetscInt    depth = 1, testNum = 0, p;
195c4762a1bSJed Brown   PetscMPIInt rank;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   PetscFunctionBegin;
1989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
199dd400576SPatrick Sanan   if (rank == 0) {
200c4762a1bSJed Brown     switch (testNum) {
2019371c9d4SSatish Balay     case 0: {
202c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
203c4762a1bSJed Brown       PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
204c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 4, 3, 5, 3, 4, 6, 5};
205c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
206c4762a1bSJed 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};
207c4762a1bSJed Brown       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
208c4762a1bSJed Brown 
2099566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
210*48a46eb9SPierre Jolivet       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
2119371c9d4SSatish Balay     } break;
2129371c9d4SSatish Balay     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
213c4762a1bSJed Brown     }
214c4762a1bSJed Brown   } else {
215c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
216c4762a1bSJed Brown 
2179566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2189566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
219c4762a1bSJed Brown   }
220c4762a1bSJed Brown   PetscFunctionReturn(0);
221c4762a1bSJed Brown }
222c4762a1bSJed Brown 
2239371c9d4SSatish Balay PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm) {
224c4762a1bSJed Brown   PetscInt    depth = 1, p;
225c4762a1bSJed Brown   PetscMPIInt rank;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   PetscFunctionBegin;
2289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
229dd400576SPatrick Sanan   if (rank == 0) {
230c4762a1bSJed Brown     switch (testNum) {
2319371c9d4SSatish Balay     case 0: {
232c4762a1bSJed Brown       PetscInt    numPoints[2]        = {6, 2};
233c4762a1bSJed Brown       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
234c4762a1bSJed Brown       PetscInt    cones[8]            = {2, 3, 4, 5, 3, 6, 7, 4};
235c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
236c4762a1bSJed 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};
237c4762a1bSJed Brown       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
240*48a46eb9SPierre Jolivet       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
2419371c9d4SSatish Balay     } break;
2429371c9d4SSatish Balay     case 1: {
243c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
244c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 3, 0, 0, 0, 0, 0};
245c4762a1bSJed Brown       PetscInt    cones[7]            = {2, 3, 4, 5, 3, 6, 4};
246c4762a1bSJed Brown       PetscInt    coneOrientations[7] = {0, 0, 0, 0, 0, 0, 0};
247c4762a1bSJed Brown       PetscScalar vertexCoords[10]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0};
248c4762a1bSJed Brown       PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
249c4762a1bSJed Brown 
2509566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
251*48a46eb9SPierre Jolivet       for (p = 0; p < 5; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
2529371c9d4SSatish Balay     } break;
2539371c9d4SSatish Balay     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
254c4762a1bSJed Brown     }
255c4762a1bSJed Brown   } else {
256c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
257c4762a1bSJed Brown 
2589566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2599566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
260c4762a1bSJed Brown   }
261c4762a1bSJed Brown   PetscFunctionReturn(0);
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
2649371c9d4SSatish Balay PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm) {
265c4762a1bSJed Brown   PetscInt    depth = 1, testNum = 0, p;
266c4762a1bSJed Brown   PetscMPIInt rank;
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   PetscFunctionBegin;
2699566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
270dd400576SPatrick Sanan   if (rank == 0) {
271c4762a1bSJed Brown     switch (testNum) {
2729371c9d4SSatish Balay     case 0: {
273c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
274c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
275c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8};
276c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
2779371c9d4SSatish Balay       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, -0.5, 0.0, 1.0, 0.0, 0.0, 1.0, 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};
278c4762a1bSJed Brown       PetscInt    markerPoints[16]     = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1};
279c4762a1bSJed Brown 
2809566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
281*48a46eb9SPierre Jolivet       for (p = 0; p < 8; ++p) PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1]));
2829371c9d4SSatish Balay     } break;
2839371c9d4SSatish Balay     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
284c4762a1bSJed Brown     }
285c4762a1bSJed Brown   } else {
286c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
287c4762a1bSJed Brown 
2889566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
2899566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(dm, "marker"));
290c4762a1bSJed Brown   }
291c4762a1bSJed Brown   PetscFunctionReturn(0);
292c4762a1bSJed Brown }
293c4762a1bSJed Brown 
2949371c9d4SSatish Balay PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm) {
295c4762a1bSJed Brown   PetscBool useGenerator = user->useGenerator;
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   PetscFunctionBegin;
2989566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2999566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
300b5a892a1SMatthew G. Knepley   if (!useGenerator) {
30130602db0SMatthew G. Knepley     PetscInt  dim     = user->dim;
30230602db0SMatthew G. Knepley     PetscBool simplex = user->simplex;
30330602db0SMatthew G. Knepley 
3049566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dm, dim));
305c4762a1bSJed Brown     switch (dim) {
306c4762a1bSJed Brown     case 2:
30730602db0SMatthew G. Knepley       if (simplex) {
3089566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_2D(comm, *dm));
309c4762a1bSJed Brown       } else {
3109566063dSJacob Faibussowitsch         PetscCall(CreateQuad_2D(comm, testNum, *dm));
311c4762a1bSJed Brown       }
312c4762a1bSJed Brown       break;
313c4762a1bSJed Brown     case 3:
31430602db0SMatthew G. Knepley       if (simplex) {
3159566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_3D(comm, *dm));
316c4762a1bSJed Brown       } else {
3179566063dSJacob Faibussowitsch         PetscCall(CreateHex_3D(comm, *dm));
318c4762a1bSJed Brown       }
319c4762a1bSJed Brown       break;
3209371c9d4SSatish Balay     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
321c4762a1bSJed Brown     }
322c4762a1bSJed Brown   }
3239566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
3249566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Mesh"));
3259566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
326c4762a1bSJed Brown   PetscFunctionReturn(0);
327c4762a1bSJed Brown }
328c4762a1bSJed Brown 
3299371c9d4SSatish Balay int main(int argc, char **argv) {
33030602db0SMatthew G. Knepley   DM     dm;
331c4762a1bSJed Brown   AppCtx user; /* user-defined work context */
332c4762a1bSJed Brown 
333327415f7SBarry Smith   PetscFunctionBeginUser;
3349566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3359566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
3369566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm));
3379566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
3389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
339b122ec5aSJacob Faibussowitsch   return 0;
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown /*TEST
343b5a892a1SMatthew G. Knepley   testset:
344b5a892a1SMatthew G. Knepley     args: -dm_plex_interpolate_pre -dm_plex_check_all
345c4762a1bSJed Brown 
346c4762a1bSJed Brown     # Two cell test meshes 0-7
347c4762a1bSJed Brown     test:
348c4762a1bSJed Brown       suffix: 0
34967d101c4SVaclav Hapla       args: -dim 2 -dm_view ascii::ascii_info_detail
350c4762a1bSJed Brown     test:
351c4762a1bSJed Brown       suffix: 1
352c4762a1bSJed Brown       nsize: 2
353e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail
354c4762a1bSJed Brown     test:
355c4762a1bSJed Brown       suffix: 2
35630602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
357c4762a1bSJed Brown     test:
358c4762a1bSJed Brown       suffix: 3
359c4762a1bSJed Brown       nsize: 2
360e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
361c4762a1bSJed Brown     test:
362c4762a1bSJed Brown       suffix: 4
363c4762a1bSJed Brown       args: -dim 3 -dm_view ascii::ascii_info_detail
364c4762a1bSJed Brown     test:
365c4762a1bSJed Brown       suffix: 5
366c4762a1bSJed Brown       nsize: 2
367e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail
368c4762a1bSJed Brown     test:
369c4762a1bSJed Brown       suffix: 6
37030602db0SMatthew G. Knepley       args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
371c4762a1bSJed Brown     test:
372c4762a1bSJed Brown       suffix: 7
373c4762a1bSJed Brown       nsize: 2
374e600fa54SMatthew G. Knepley       args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
375c4762a1bSJed Brown     # 2D Hybrid Mesh 8
376c4762a1bSJed Brown     test:
377c4762a1bSJed Brown       suffix: 8
37830602db0SMatthew G. Knepley       args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail
379b5a892a1SMatthew G. Knepley 
380b5a892a1SMatthew G. Knepley   testset:
381b5a892a1SMatthew G. Knepley     args: -dm_plex_check_all
382b5a892a1SMatthew G. Knepley 
383da9060c4SMatthew G. Knepley     # Reference cells
384da9060c4SMatthew G. Knepley     test:
385da9060c4SMatthew G. Knepley       suffix: 12
386b5a892a1SMatthew G. Knepley       args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all
387c4762a1bSJed Brown     # TetGen meshes 9-10
388c4762a1bSJed Brown     test:
389c4762a1bSJed Brown       suffix: 9
390c4762a1bSJed Brown       requires: triangle
39130602db0SMatthew G. Knepley       args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0
392c4762a1bSJed Brown     test:
393c4762a1bSJed Brown       suffix: 10
394c4762a1bSJed Brown       requires: ctetgen
39530602db0SMatthew G. Knepley       args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0
396c4762a1bSJed Brown     # Cubit meshes 11
397c4762a1bSJed Brown     test:
398c4762a1bSJed Brown       suffix: 11
399c4762a1bSJed Brown       requires: exodusii
40030602db0SMatthew 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
401c4762a1bSJed Brown 
402c4762a1bSJed Brown TEST*/
403