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