1cf4936a9SMatthew G. Knepley static char help[] = "PETSc Annual Meeting 2025: Meshing Tutorial.\n\n\n";
2cf4936a9SMatthew G. Knepley
3cf4936a9SMatthew G. Knepley #include <petsc.h>
4cf4936a9SMatthew G. Knepley
CreateMesh(MPI_Comm comm,DM * dm)5cf4936a9SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
6cf4936a9SMatthew G. Knepley {
7cf4936a9SMatthew G. Knepley PetscFunctionBeginUser;
8cf4936a9SMatthew G. Knepley PetscCall(DMCreate(comm, dm));
9cf4936a9SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX));
10cf4936a9SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm));
11cf4936a9SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
12cf4936a9SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS);
13cf4936a9SMatthew G. Knepley }
14cf4936a9SMatthew G. Knepley
main(int argc,char ** argv)15cf4936a9SMatthew G. Knepley int main(int argc, char **argv)
16cf4936a9SMatthew G. Knepley {
17cf4936a9SMatthew G. Knepley DM dm;
18cf4936a9SMatthew G. Knepley
19cf4936a9SMatthew G. Knepley PetscFunctionBeginUser;
20cf4936a9SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help));
21cf4936a9SMatthew G. Knepley PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
22cf4936a9SMatthew G. Knepley PetscCall(DMDestroy(&dm));
23cf4936a9SMatthew G. Knepley PetscCall(PetscFinalize());
24cf4936a9SMatthew G. Knepley return PETSC_SUCCESS;
25cf4936a9SMatthew G. Knepley }
26cf4936a9SMatthew G. Knepley
27cf4936a9SMatthew G. Knepley /*TEST
28cf4936a9SMatthew G. Knepley
29cf4936a9SMatthew G. Knepley # Draw a square with X, use with -draw_pause -1
30cf4936a9SMatthew G. Knepley test:
31cf4936a9SMatthew G. Knepley suffix: 0
32cf4936a9SMatthew G. Knepley requires: triangle x
33cf4936a9SMatthew G. Knepley args: -dm_view draw -draw_pause 0
34*3886731fSPierre Jolivet output_file: output/empty.out
35cf4936a9SMatthew G. Knepley
36cf4936a9SMatthew G. Knepley # Draw a square with PyVista
37cf4936a9SMatthew G. Knepley test:
38cf4936a9SMatthew G. Knepley suffix: 1
39cf4936a9SMatthew G. Knepley requires: triangle pyvista
40cf4936a9SMatthew G. Knepley args: -dm_view pyvista
41*3886731fSPierre Jolivet output_file: output/empty.out
42cf4936a9SMatthew G. Knepley
43cf4936a9SMatthew G. Knepley # Refine the square
44cf4936a9SMatthew G. Knepley test:
45cf4936a9SMatthew G. Knepley suffix: 2
46cf4936a9SMatthew G. Knepley requires: triangle pyvista
47cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_refine 1
48*3886731fSPierre Jolivet output_file: output/empty.out
49cf4936a9SMatthew G. Knepley
50cf4936a9SMatthew G. Knepley # Refine the square three times
51cf4936a9SMatthew G. Knepley test:
52cf4936a9SMatthew G. Knepley suffix: 3
53cf4936a9SMatthew G. Knepley requires: triangle pyvista
54cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_refine 3
55*3886731fSPierre Jolivet output_file: output/empty.out
56cf4936a9SMatthew G. Knepley
57cf4936a9SMatthew G. Knepley # Refine the cube three times
58cf4936a9SMatthew G. Knepley test:
59cf4936a9SMatthew G. Knepley suffix: 4
60cf4936a9SMatthew G. Knepley requires: ctetgen pyvista
61cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_dim 3 -dm_refine 3
62*3886731fSPierre Jolivet output_file: output/empty.out
63cf4936a9SMatthew G. Knepley
64cf4936a9SMatthew G. Knepley # Draw a sphere with PyVista (all we get is an icosahedron)
65cf4936a9SMatthew G. Knepley test:
66cf4936a9SMatthew G. Knepley suffix: 5
67cf4936a9SMatthew G. Knepley requires: pyvista
68cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere
69*3886731fSPierre Jolivet output_file: output/empty.out
70cf4936a9SMatthew G. Knepley
71cf4936a9SMatthew G. Knepley # Refine the sphere three times
72cf4936a9SMatthew G. Knepley test:
73cf4936a9SMatthew G. Knepley suffix: 6
74cf4936a9SMatthew G. Knepley requires: pyvista
75cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere -dm_refine 3
76*3886731fSPierre Jolivet output_file: output/empty.out
77cf4936a9SMatthew G. Knepley
78cf4936a9SMatthew G. Knepley # Show the 3-sphere
79cf4936a9SMatthew G. Knepley test:
80cf4936a9SMatthew G. Knepley suffix: 7
81cf4936a9SMatthew G. Knepley args: -dm_view -dm_plex_shape sphere -dm_plex_dim 3
82cf4936a9SMatthew G. Knepley
83cf4936a9SMatthew G. Knepley # Extrude the sphere
84cf4936a9SMatthew G. Knepley test:
85cf4936a9SMatthew G. Knepley suffix: 8
86cf4936a9SMatthew G. Knepley requires: pyvista
87cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere -dm_refine 1 \
88cf4936a9SMatthew G. Knepley -dm_plex_transform_type extrude -dm_plex_transform_extrude_layers 3 \
89cf4936a9SMatthew G. Knepley -dm_plex_transform_extrude_use_tensor 0 -dm_plex_transform_extrude_thickness 0.3
90*3886731fSPierre Jolivet output_file: output/empty.out
91cf4936a9SMatthew G. Knepley
92cf4936a9SMatthew G. Knepley # Extrude the sphere with cutaway
93cf4936a9SMatthew G. Knepley test:
94cf4936a9SMatthew G. Knepley suffix: 9
95cf4936a9SMatthew G. Knepley requires: pyvista
96cf4936a9SMatthew G. Knepley args: -dm_view pyvista -view_pyvista_clip 1.,0.,0. -dm_plex_shape sphere -dm_refine 1 \
97cf4936a9SMatthew G. Knepley -dm_plex_transform_type extrude -dm_plex_transform_extrude_layers 3 \
98cf4936a9SMatthew G. Knepley -dm_plex_transform_extrude_use_tensor 0 -dm_plex_transform_extrude_thickness 0.3
99*3886731fSPierre Jolivet output_file: output/empty.out
100cf4936a9SMatthew G. Knepley
101cf4936a9SMatthew G. Knepley # Extrude the refined sphere
102cf4936a9SMatthew G. Knepley test:
103cf4936a9SMatthew G. Knepley suffix: 10
104cf4936a9SMatthew G. Knepley requires: pyvista
105cf4936a9SMatthew G. Knepley args: -dm_view pyvista -view_pyvista_clip 1.,0.,0. -dm_plex_shape sphere -dm_refine 3 \
106cf4936a9SMatthew G. Knepley -dm_plex_option_phases ext_ \
107cf4936a9SMatthew G. Knepley -ext_dm_refine 1 -ext_dm_plex_transform_type extrude -ext_dm_plex_transform_extrude_layers 3 \
108cf4936a9SMatthew G. Knepley -ext_dm_plex_transform_extrude_use_tensor 0 -ext_dm_plex_transform_extrude_thickness 0.5
109*3886731fSPierre Jolivet output_file: output/empty.out
110cf4936a9SMatthew G. Knepley
111cf4936a9SMatthew G. Knepley # Extrude the refined sphere
112cf4936a9SMatthew G. Knepley test:
113cf4936a9SMatthew G. Knepley suffix: 11
114cf4936a9SMatthew G. Knepley requires: pyvista
115cf4936a9SMatthew G. Knepley args: -dm_view pyvista -view_pyvista_clip 1.,0.,0. -dm_plex_shape sphere -dm_refine 3 \
116cf4936a9SMatthew G. Knepley -dm_plex_option_phases ext_,ref_ \
117cf4936a9SMatthew G. Knepley -ext_dm_refine 1 -ext_dm_plex_transform_type extrude -ext_dm_plex_transform_extrude_layers 3 \
118cf4936a9SMatthew G. Knepley -ext_dm_plex_transform_extrude_use_tensor 0 -ext_dm_plex_transform_extrude_thickness 0.5 \
119cf4936a9SMatthew G. Knepley -ref_dm_refine 1
120*3886731fSPierre Jolivet output_file: output/empty.out
121cf4936a9SMatthew G. Knepley
122cf4936a9SMatthew G. Knepley # Refine the Schwartz surface
123cf4936a9SMatthew G. Knepley test:
124cf4936a9SMatthew G. Knepley suffix: 12
125cf4936a9SMatthew G. Knepley requires: pyvista
126cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape schwarz_p -dm_plex_tps_extent 3,2 -dm_plex_tps_refine 2
127*3886731fSPierre Jolivet output_file: output/empty.out
128cf4936a9SMatthew G. Knepley
129cf4936a9SMatthew G. Knepley # Refine and extrude the Schwartz surface with given thickness
130cf4936a9SMatthew G. Knepley test:
131cf4936a9SMatthew G. Knepley suffix: 13
132cf4936a9SMatthew G. Knepley requires: pyvista
133cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape schwarz_p -dm_plex_tps_extent 3,2,2 -dm_plex_tps_refine 3 \
134cf4936a9SMatthew G. Knepley -dm_plex_tps_thickness 0.4 -dm_plex_tps_layers 3
135*3886731fSPierre Jolivet output_file: output/empty.out
136cf4936a9SMatthew G. Knepley
137cf4936a9SMatthew G. Knepley # Filter the square
138cf4936a9SMatthew G. Knepley test:
139cf4936a9SMatthew G. Knepley suffix: 14
140cf4936a9SMatthew G. Knepley requires: triangle pyvista
141cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_refine 1 \
142cf4936a9SMatthew G. Knepley -dm_plex_transform_type transform_filter -dm_plex_transform_active filter_cells -dm_plex_label_filter_cells 0,1,2
143*3886731fSPierre Jolivet output_file: output/empty.out
144cf4936a9SMatthew G. Knepley
145cf4936a9SMatthew G. Knepley # Filter a cap on the sphere
146cf4936a9SMatthew G. Knepley test:
147cf4936a9SMatthew G. Knepley suffix: 15
148cf4936a9SMatthew G. Knepley requires: pyvista
149cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere -dm_refine 4 \
150cf4936a9SMatthew G. Knepley -dm_plex_option_phases filt_ \
151cf4936a9SMatthew G. Knepley -filt_dm_refine 1 -filt_dm_plex_transform_type transform_filter \
152cf4936a9SMatthew G. Knepley -filt_dm_plex_transform_active filter_cells -dm_plex_label_filter_cells 0
153*3886731fSPierre Jolivet output_file: output/empty.out
154cf4936a9SMatthew G. Knepley
155cf4936a9SMatthew G. Knepley # Filter the sphere minus a cap
156cf4936a9SMatthew G. Knepley test:
157cf4936a9SMatthew G. Knepley suffix: 16
158cf4936a9SMatthew G. Knepley requires: pyvista
159cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere -dm_refine 4 \
160cf4936a9SMatthew G. Knepley -dm_plex_option_phases filt_ \
161cf4936a9SMatthew G. Knepley -filt_dm_refine 1 -filt_dm_plex_transform_type transform_filter \
162cf4936a9SMatthew G. Knepley -filt_dm_plex_transform_active filter_cells \
163cf4936a9SMatthew G. Knepley -dm_plex_label_filter_cells 0,1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16,17,18
164*3886731fSPierre Jolivet output_file: output/empty.out
165cf4936a9SMatthew G. Knepley
166cf4936a9SMatthew G. Knepley # Convert sphere to quads
167cf4936a9SMatthew G. Knepley test:
168cf4936a9SMatthew G. Knepley suffix: 17
169cf4936a9SMatthew G. Knepley requires: pyvista
170cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_shape sphere -dm_refine 3 \
171cf4936a9SMatthew G. Knepley -dm_plex_option_phases conv_ -conv_dm_refine 1 -conv_dm_plex_transform_type refine_tobox
172*3886731fSPierre Jolivet output_file: output/empty.out
173cf4936a9SMatthew G. Knepley
174cf4936a9SMatthew G. Knepley # Load and refine the nozzle
175cf4936a9SMatthew G. Knepley test:
176cf4936a9SMatthew G. Knepley suffix: 18
177cf4936a9SMatthew G. Knepley requires: pyvista egads datafilespath
178cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_filename /Users/knepley/PETSc4/petsc/petsc-dev/share/petsc/datafiles/meshes/nozzle.stp -dm_refine 3
179*3886731fSPierre Jolivet output_file: output/empty.out
180cf4936a9SMatthew G. Knepley
181cf4936a9SMatthew G. Knepley # Load and refine the nozzle, and convert to quads
182cf4936a9SMatthew G. Knepley test:
183cf4936a9SMatthew G. Knepley suffix: 19
184cf4936a9SMatthew G. Knepley requires: pyvista egads datafilespath
185cf4936a9SMatthew G. Knepley args: -dm_view pyvista -dm_plex_filename /Users/knepley/PETSc4/petsc/petsc-dev/share/petsc/datafiles/meshes/nozzle.stp -dm_refine 3 \
186cf4936a9SMatthew G. Knepley -dm_plex_option_phases conv_ -conv_dm_refine 1 -conv_dm_plex_transform_type refine_tobox
187*3886731fSPierre Jolivet output_file: output/empty.out
188cf4936a9SMatthew G. Knepley
189cf4936a9SMatthew G. Knepley TEST*/
190