xref: /petsc/src/dm/impls/plex/tests/ex40.c (revision a2fddd78f770fa4fc19a8af67e65be331f27d92b)
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_tri
35     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_refine 2 -dm_plex_check_all
36 
37   test:
38     suffix: box_tri
39     requires: triangle
40     nsize: {{1 3 5}}
41     args: -dm_distribute -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_check_all
42 
43   test:
44     suffix: ref_quad
45     args: -dm_plex_reference_cell_domain -dm_plex_cell quadrilateral -dm_refine 2 -dm_plex_check_all
46 
47   test:
48     suffix: box_quad
49     nsize: {{1 3 5}}
50     args: -dm_distribute -dm_plex_box_faces 3,3 -dm_plex_simplex 0 -dm_refine 2 -dm_plex_check_all
51 
52   test:
53     suffix: ref_tet
54     args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_refine 2 -dm_plex_check_all
55 
56   test:
57     suffix: box_tet
58     requires: ctetgen
59     nsize: {{1 3 5}}
60     args: -dm_distribute -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 -dm_plex_check_all
61 
62   test:
63     suffix: ref_hex
64     args: -dm_plex_reference_cell_domain -dm_plex_cell hexahedron -dm_refine 2 -dm_plex_check_all
65 
66   test:
67     suffix: box_hex
68     nsize: {{1 3 5}}
69     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
70 
71   test:
72     suffix: ref_trip
73     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_refine 2 -dm_plex_check_all
74 
75   test:
76     suffix: ref_tquad
77     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quad -dm_refine 2 -dm_plex_check_all
78 
79   test:
80     suffix: ref_ttrip
81     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_refine 2 -dm_plex_check_all
82 
83   test:
84     suffix: ref_tquadp
85     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_quadrilateral_prism -dm_refine 2 -dm_plex_check_all
86 
87   test:
88     suffix: ref_tri_tobox
89     args: -dm_plex_reference_cell_domain -dm_plex_cell triangle -dm_plex_cell_refiner tobox -dm_coord_space 0 -dm_refine 2 -dm_plex_check_all
90 
91   test:
92     suffix: box_tri_tobox
93     requires: triangle
94     nsize: {{1 3 5}}
95     args: -dm_distribute -dm_plex_box_faces 3,3 -dm_plex_cell_refiner tobox -dm_coord_space 0 -dm_refine 2 -dm_plex_check_all
96 
97   test:
98     suffix: ref_tet_tobox
99     args: -dm_plex_reference_cell_domain -dm_plex_cell tetrahedron -dm_plex_cell_refiner tobox -dm_coord_space 0 -dm_refine 2 -dm_plex_check_all
100 
101   test:
102     suffix: box_tet_tobox
103     requires: ctetgen
104     nsize: {{1 3 5}}
105     args: -dm_distribute -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_cell_refiner tobox -dm_coord_space 0 -dm_refine 2 -dm_plex_check_all
106 
107   test:
108     suffix: ref_trip_tobox
109     args: -dm_plex_reference_cell_domain -dm_plex_cell triangular_prism -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
110 
111   test:
112     suffix: ref_ttrip_tobox
113     args: -dm_plex_reference_cell_domain -dm_plex_cell tensor_triangular_prism -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
114 
115 TEST*/
116