xref: /petsc/src/dm/impls/plex/tests/ex40.c (revision daad07d386296cdcbb87925ef5f1432ee4a24ec4)
1 static const char help[] = "Tests for Plex transforms, including regular refinement";
2 
3 #include <petscdmplex.h>
4 #include <petscsf.h>
5 
6 #include <petsc/private/dmpleximpl.h>
7 
8 static PetscErrorCode LabelPoints(DM dm)
9 {
10   DMLabel        label;
11   PetscInt       pStart, pEnd, p;
12   PetscBool      flg = PETSC_FALSE;
13 
14   PetscFunctionBegin;
15   PetscCall(PetscOptionsGetBool(NULL, NULL, "-label_mesh", &flg, NULL));
16   if (!flg) PetscFunctionReturn(0);
17   PetscCall(DMCreateLabel(dm, "test"));
18   PetscCall(DMGetLabel(dm, "test", &label));
19   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
20   for (p = pStart; p < pEnd; ++p) {
21     PetscCall(DMLabelSetValue(label, p, p));
22   }
23   PetscFunctionReturn(0);
24 }
25 
26 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
27 {
28   PetscFunctionBegin;
29   PetscCall(DMCreate(comm, dm));
30   PetscCall(DMSetType(*dm, DMPLEX));
31   PetscCall(DMSetFromOptions(*dm));
32   PetscCall(LabelPoints(*dm));
33   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, "post_label_"));
34   PetscCall(DMSetFromOptions(*dm));
35   PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *dm, NULL));
36   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
37   PetscFunctionReturn(0);
38 }
39 
40 int main(int argc, char **argv)
41 {
42   DM             dm;
43 
44   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
45   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
46   PetscCall(DMDestroy(&dm));
47   PetscCall(PetscFinalize());
48   return 0;
49 }
50 
51 /*TEST
52   test:
53     suffix: ref_seg
54     args: -dm_plex_reference_cell_domain -dm_plex_cell segment -dm_refine 1 -dm_plex_check_all
55 
56   test:
57     suffix: ref_tri
58     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2 -dm_plex_check_all
59 
60   test:
61     suffix: box_tri
62     requires: triangle
63     nsize: {{1 3 5}}
64     args: -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_check_all
65 
66   test:
67     suffix: ref_quad
68     args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -dm_refine 2 -dm_plex_check_all
69 
70   test:
71     suffix: box_quad
72     nsize: {{1 3 5}}
73     args: -dm_plex_box_faces 3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all
74 
75   test:
76     suffix: ref_tet
77     args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_refine 2 -dm_plex_check_all
78 
79   test:
80     suffix: box_tet
81     requires: ctetgen
82     nsize: {{1 3 5}}
83     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 -dm_plex_check_all
84 
85   test:
86     suffix: ref_hex
87     args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -dm_refine 2 -dm_plex_check_all
88 
89   test:
90     suffix: box_hex
91     nsize: {{1 3 5}}
92     args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all
93 
94   test:
95     suffix: ref_trip
96     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_refine 2 -dm_plex_check_all
97 
98   test:
99     suffix: ref_tquad
100     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -dm_refine 2 -dm_plex_check_all
101 
102   test:
103     suffix: ref_ttrip
104     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_refine 2 -dm_plex_check_all
105 
106   test:
107     suffix: ref_tquadp
108     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -dm_refine 2 -dm_plex_check_all
109 
110   test:
111     suffix: ref_pyramid
112     args: -dm_plex_reference_cell_domain -dm_plex_cell pyramid -dm_refine 2 -dm_plex_check_all
113 
114   testset:
115     args: -dm_coord_space 0 -dm_plex_transform_type refine_tobox -dm_plex_check_all
116 
117     test:
118       suffix: ref_tri_tobox
119       args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2
120 
121     test:
122       suffix: box_tri_tobox
123       requires: triangle
124       nsize: {{1 3 5}}
125       args: -dm_plex_box_faces 3,3 -dm_refine 2
126 
127     test:
128       suffix: ref_tet_tobox
129       args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_refine 2
130 
131     test:
132       suffix: box_tet_tobox
133       requires: ctetgen
134       nsize: {{1 3 5}}
135       args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2
136 
137     test:
138       suffix: ref_trip_tobox
139       args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_refine 2
140 
141     test:
142       suffix: ref_ttrip_tobox
143       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_refine 2
144 
145     test:
146       suffix: ref_tquadp_tobox
147       args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -dm_refine 2
148 
149   testset:
150     args: -dm_coord_space 0 -label_mesh -post_label_dm_extrude 2 -post_label_dm_plex_check_all -dm_view ::ascii_info_detail
151 
152     test:
153       suffix: extrude_quad
154       args: -dm_plex_simplex 0
155 
156 TEST*/
157