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