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