xref: /petsc/src/dm/impls/plex/tests/ex40.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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) PetscCall(DMLabelSetValue(label, p, p));
21   PetscFunctionReturn(0);
22 }
23 
24 static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
25 {
26   PetscFunctionBegin;
27   PetscCall(DMCreate(comm, dm));
28   PetscCall(DMSetType(*dm, DMPLEX));
29   PetscCall(DMSetFromOptions(*dm));
30   PetscCall(LabelPoints(*dm));
31   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "post_label_"));
32   PetscCall(DMSetFromOptions(*dm));
33   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, NULL));
34   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
35   PetscFunctionReturn(0);
36 }
37 
38 int main(int argc, char **argv)
39 {
40   DM dm;
41 
42   PetscFunctionBeginUser;
43   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
44   PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
45   PetscCall(DMDestroy(&dm));
46   PetscCall(PetscFinalize());
47   return 0;
48 }
49 
50 /*TEST
51   test:
52     suffix: ref_seg
53     args: -dm_plex_reference_cell_domain -dm_plex_cell segment -dm_refine 1 -dm_plex_check_all
54 
55   test:
56     suffix: ref_tri
57     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2 -dm_plex_check_all
58 
59   test:
60     suffix: box_tri
61     requires: triangle parmetis
62     nsize: {{1 3 5}}
63     args: -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_check_all -petscpartitioner_type parmetis
64 
65   test:
66     suffix: ref_quad
67     args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -dm_refine 2 -dm_plex_check_all
68 
69   test:
70     suffix: box_quad
71     nsize: {{1 3 5}}
72     requires: parmetis
73     args: -dm_plex_box_faces 3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all -petscpartitioner_type parmetis
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 parmetis
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 -petscpartitioner_type parmetis
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     requires: parmetis
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 -petscpartitioner_type parmetis
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 parmetis
125       nsize: {{1 3 5}}
126       args: -dm_plex_box_faces 3,3 -dm_refine 2 -petscpartitioner_type parmetis
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 parmetis
135       nsize: {{1 3 5}}
136       args: -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 -petscpartitioner_type parmetis
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