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