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