xref: /petsc/src/dm/impls/plex/tests/ex7.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 static char help[] = "Tests for mesh interpolation\n\n";
2 
3 #include <petscdmplex.h>
4 
5 /* List of test meshes
6 
7 Triangle
8 --------
9 Test 0:
10 Two triangles sharing a face
11 
12         4
13       / | \
14      /  |  \
15     /   |   \
16    2  0 | 1  5
17     \   |   /
18      \  |  /
19       \ | /
20         3
21 
22 should become
23 
24         4
25       / | \
26      8  |  9
27     /   |   \
28    2  0 7 1  5
29     \   |   /
30      6  |  10
31       \ | /
32         3
33 
34 Tetrahedron
35 -----------
36 Test 0:
37 Two tets sharing a face
38 
39  cell   5 _______    cell
40  0    / | \      \       1
41      /  |  \      \
42     /   |   \      \
43    2----|----4-----6
44     \   |   /      /
45      \  |  /     /
46       \ | /      /
47         3-------
48 
49 should become
50 
51  cell   5 _______    cell
52  0    / | \      \       1
53      /  |  \      \
54    17   |   18 13  22
55    / 8 19 10 \      \
56   /     |     \      \
57  2---14-|------4--20--6
58   \     |     /      /
59    \ 9  | 7  /      /
60    16   |   15 11  21
61      \  |  /      /
62       \ | /      /
63         3-------
64 
65 Quadrilateral
66 -------------
67 Test 0:
68 Two quads sharing a face
69 
70    5-------4-------7
71    |       |       |
72    |   0   |   1   |
73    |       |       |
74    2-------3-------6
75 
76 should become
77 
78    5--10---4--14---7
79    |       |       |
80   11   0   9   1  13
81    |       |       |
82    2---8---3--12---6
83 
84 Test 1:
85 A quad and a triangle sharing a face
86 
87    5-------4
88    |       | \
89    |   0   |  \
90    |       | 1 \
91    2-------3----6
92 
93 should become
94 
95    5---9---4
96    |       | \
97   10   0   8  12
98    |       | 1 \
99    2---7---3-11-6
100 
101 Hexahedron
102 ----------
103 Test 0:
104 Two hexes sharing a face
105 
106 cell   9-------------8-------------13 cell
107 0     /|            /|            /|     1
108      / |   15      / |   21      / |
109     /  |          /  |          /  |
110    6-------------7-------------12  |
111    |   |     18  |   |     24  |   |
112    |   |         |   |         |   |
113    |19 |         |17 |         |23 |
114    |   |  16     |   |   22    |   |
115    |   5---------|---4---------|---11
116    |  /          |  /          |  /
117    | /   14      | /    20     | /
118    |/            |/            |/
119    2-------------3-------------10
120 
121 should become
122 
123 cell   9-----31------8-----42------13 cell
124 0     /|            /|            /|     1
125     32 |   15      30|   21      41|
126     /  |          /  |          /  |
127    6-----29------7-----40------12  |
128    |   |     17  |   |     23  |   |
129    |  35         |  36         |   44
130    |19 |         |18 |         |24 |
131   34   |  16    33   |   22   43   |
132    |   5-----26--|---4-----37--|---11
133    |  /          |  /          |  /
134    | 25   14     | 27    20    | 38
135    |/            |/            |/
136    2-----28------3-----39------10
137 
138 */
139 
140 typedef struct {
141   PetscInt  testNum;      /* Indicates the mesh to create */
142   PetscInt  dim;          /* The topological mesh dimension */
143   PetscBool simplex;      /* Use simplices or hexes */
144   PetscBool useGenerator; /* Construct mesh with a mesh generator */
145 } AppCtx;
146 
147 PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) {
148   PetscFunctionBegin;
149   options->testNum      = 0;
150   options->dim          = 2;
151   options->simplex      = PETSC_TRUE;
152   options->useGenerator = PETSC_FALSE;
153 
154   PetscOptionsBegin(comm, "", "Meshing Interpolation Test Options", "DMPLEX");
155   PetscCall(PetscOptionsBoundedInt("-testnum", "The mesh to create", "ex7.c", options->testNum, &options->testNum, NULL, 0));
156   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex7.c", options->dim, &options->dim, NULL, 1, 3));
157   PetscCall(PetscOptionsBool("-simplex", "Use simplices if true, otherwise hexes", "ex7.c", options->simplex, &options->simplex, NULL));
158   PetscCall(PetscOptionsBool("-use_generator", "Use a mesh generator to build the mesh", "ex7.c", options->useGenerator, &options->useGenerator, NULL));
159   PetscOptionsEnd();
160   PetscFunctionReturn(0);
161 }
162 
163 PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM dm) {
164   PetscInt    depth = 1, testNum = 0, p;
165   PetscMPIInt rank;
166 
167   PetscFunctionBegin;
168   PetscCallMPI(MPI_Comm_rank(comm, &rank));
169   if (rank == 0) {
170     switch (testNum) {
171     case 0: {
172       PetscInt    numPoints[2]        = {4, 2};
173       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
174       PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
175       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
176       PetscScalar vertexCoords[8]     = {-0.5, 0.5, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5};
177       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
178 
179       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
180       for (p = 0; p < 4; ++p) { PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); }
181     } break;
182     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
183     }
184   } else {
185     PetscInt numPoints[2] = {0, 0};
186 
187     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
188     PetscCall(DMCreateLabel(dm, "marker"));
189   }
190   PetscFunctionReturn(0);
191 }
192 
193 PetscErrorCode CreateSimplex_3D(MPI_Comm comm, DM dm) {
194   PetscInt    depth = 1, testNum = 0, p;
195   PetscMPIInt rank;
196 
197   PetscFunctionBegin;
198   PetscCallMPI(MPI_Comm_rank(comm, &rank));
199   if (rank == 0) {
200     switch (testNum) {
201     case 0: {
202       PetscInt    numPoints[2]        = {5, 2};
203       PetscInt    coneSize[23]        = {4, 4, 0, 0, 0, 0, 0};
204       PetscInt    cones[8]            = {2, 4, 3, 5, 3, 4, 6, 5};
205       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
206       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};
207       PetscInt    markerPoints[8]     = {2, 1, 3, 1, 4, 1, 5, 1};
208 
209       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
210       for (p = 0; p < 4; ++p) { PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); }
211     } break;
212     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
213     }
214   } else {
215     PetscInt numPoints[2] = {0, 0};
216 
217     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
218     PetscCall(DMCreateLabel(dm, "marker"));
219   }
220   PetscFunctionReturn(0);
221 }
222 
223 PetscErrorCode CreateQuad_2D(MPI_Comm comm, PetscInt testNum, DM dm) {
224   PetscInt    depth = 1, p;
225   PetscMPIInt rank;
226 
227   PetscFunctionBegin;
228   PetscCallMPI(MPI_Comm_rank(comm, &rank));
229   if (rank == 0) {
230     switch (testNum) {
231     case 0: {
232       PetscInt    numPoints[2]        = {6, 2};
233       PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
234       PetscInt    cones[8]            = {2, 3, 4, 5, 3, 6, 7, 4};
235       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
236       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};
237       PetscInt    markerPoints[12]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1};
238 
239       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
240       for (p = 0; p < 6; ++p) { PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); }
241     } break;
242     case 1: {
243       PetscInt    numPoints[2]        = {5, 2};
244       PetscInt    coneSize[7]         = {4, 3, 0, 0, 0, 0, 0};
245       PetscInt    cones[7]            = {2, 3, 4, 5, 3, 6, 4};
246       PetscInt    coneOrientations[7] = {0, 0, 0, 0, 0, 0, 0};
247       PetscScalar vertexCoords[10]    = {-0.5, 0.0, 0.0, 0.0, 0.0, 1.0, -0.5, 1.0, 0.5, 0.0};
248       PetscInt    markerPoints[10]    = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1};
249 
250       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
251       for (p = 0; p < 5; ++p) { PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); }
252     } break;
253     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
254     }
255   } else {
256     PetscInt numPoints[2] = {0, 0};
257 
258     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
259     PetscCall(DMCreateLabel(dm, "marker"));
260   }
261   PetscFunctionReturn(0);
262 }
263 
264 PetscErrorCode CreateHex_3D(MPI_Comm comm, DM dm) {
265   PetscInt    depth = 1, testNum = 0, p;
266   PetscMPIInt rank;
267 
268   PetscFunctionBegin;
269   PetscCallMPI(MPI_Comm_rank(comm, &rank));
270   if (rank == 0) {
271     switch (testNum) {
272     case 0: {
273       PetscInt    numPoints[2]         = {12, 2};
274       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
275       PetscInt    cones[16]            = {2, 5, 4, 3, 6, 7, 8, 9, 3, 4, 11, 10, 7, 12, 13, 8};
276       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
277       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};
278       PetscInt    markerPoints[16]     = {2, 1, 3, 1, 4, 1, 5, 1, 6, 1, 7, 1, 8, 1, 9, 1};
279 
280       PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
281       for (p = 0; p < 8; ++p) { PetscCall(DMSetLabelValue(dm, "marker", markerPoints[p * 2], markerPoints[p * 2 + 1])); }
282     } break;
283     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
284     }
285   } else {
286     PetscInt numPoints[2] = {0, 0};
287 
288     PetscCall(DMPlexCreateFromDAG(dm, depth, numPoints, NULL, NULL, NULL, NULL));
289     PetscCall(DMCreateLabel(dm, "marker"));
290   }
291   PetscFunctionReturn(0);
292 }
293 
294 PetscErrorCode CreateMesh(MPI_Comm comm, PetscInt testNum, AppCtx *user, DM *dm) {
295   PetscBool useGenerator = user->useGenerator;
296 
297   PetscFunctionBegin;
298   PetscCall(DMCreate(comm, dm));
299   PetscCall(DMSetType(*dm, DMPLEX));
300   if (!useGenerator) {
301     PetscInt  dim     = user->dim;
302     PetscBool simplex = user->simplex;
303 
304     PetscCall(DMSetDimension(*dm, dim));
305     switch (dim) {
306     case 2:
307       if (simplex) {
308         PetscCall(CreateSimplex_2D(comm, *dm));
309       } else {
310         PetscCall(CreateQuad_2D(comm, testNum, *dm));
311       }
312       break;
313     case 3:
314       if (simplex) {
315         PetscCall(CreateSimplex_3D(comm, *dm));
316       } else {
317         PetscCall(CreateHex_3D(comm, *dm));
318       }
319       break;
320     default: SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
321     }
322   }
323   PetscCall(DMSetFromOptions(*dm));
324   PetscCall(PetscObjectSetName((PetscObject)*dm, "Interpolated Mesh"));
325   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
326   PetscFunctionReturn(0);
327 }
328 
329 int main(int argc, char **argv) {
330   DM     dm;
331   AppCtx user; /* user-defined work context */
332 
333   PetscFunctionBeginUser;
334   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
335   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
336   PetscCall(CreateMesh(PETSC_COMM_WORLD, user.testNum, &user, &dm));
337   PetscCall(DMDestroy(&dm));
338   PetscCall(PetscFinalize());
339   return 0;
340 }
341 
342 /*TEST
343   testset:
344     args: -dm_plex_interpolate_pre -dm_plex_check_all
345 
346     # Two cell test meshes 0-7
347     test:
348       suffix: 0
349       args: -dim 2 -dm_view ascii::ascii_info_detail
350     test:
351       suffix: 1
352       nsize: 2
353       args: -petscpartitioner_type simple -dim 2 -dm_view ascii::ascii_info_detail
354     test:
355       suffix: 2
356       args: -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
357     test:
358       suffix: 3
359       nsize: 2
360       args: -petscpartitioner_type simple -dim 2 -simplex 0 -dm_view ascii::ascii_info_detail
361     test:
362       suffix: 4
363       args: -dim 3 -dm_view ascii::ascii_info_detail
364     test:
365       suffix: 5
366       nsize: 2
367       args: -petscpartitioner_type simple -dim 3 -dm_view ascii::ascii_info_detail
368     test:
369       suffix: 6
370       args: -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
371     test:
372       suffix: 7
373       nsize: 2
374       args: -petscpartitioner_type simple -dim 3 -simplex 0 -dm_view ascii::ascii_info_detail
375     # 2D Hybrid Mesh 8
376     test:
377       suffix: 8
378       args: -dim 2 -simplex 0 -testnum 1 -dm_view ascii::ascii_info_detail
379 
380   testset:
381     args: -dm_plex_check_all
382 
383     # Reference cells
384     test:
385       suffix: 12
386       args: -use_generator -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_plex_check_all
387     # TetGen meshes 9-10
388     test:
389       suffix: 9
390       requires: triangle
391       args: -use_generator -dm_view ascii::ascii_info_detail -dm_coord_space 0
392     test:
393       suffix: 10
394       requires: ctetgen
395       args: -use_generator -dm_plex_dim 3 -dm_view ascii::ascii_info_detail -dm_coord_space 0
396     # Cubit meshes 11
397     test:
398       suffix: 11
399       requires: exodusii
400       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
401 
402 TEST*/
403