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