xref: /petsc/src/dm/impls/plex/tests/ex40.c (revision 2da392cc7c10228af19ad9843ce5155178acb644)
1 static const char help[] = "Tests for regular refinement";
2 
3 /* TODO
4   - Add in simplex-to-hex tests
5 */
6 
7 #include <petscdmplex.h>
8 #include <petscsf.h>
9 
10 #include <petsc/private/dmpleximpl.h>
11 
12 typedef struct {
13   DMPolytopeType refCell; /* Use the reference cell */
14   PetscInt       dim;     /* The topological dimension */
15   PetscBool      simplex; /* Flag for simplices */
16 } AppCtx;
17 
18 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
19 {
20   PetscErrorCode ierr;
21 
22   PetscFunctionBeginUser;
23   options->refCell = DM_POLYTOPE_UNKNOWN;
24   options->dim     = 2;
25   options->simplex = PETSC_TRUE;
26 
27   ierr = PetscOptionsBegin(comm, "", "Parallel Mesh Orientation Options", "DMPLEX");CHKERRQ(ierr);
28   ierr = PetscOptionsInt("-dim", "The topological dimension", "ex40.c", options->dim, &options->dim, NULL);CHKERRQ(ierr);
29   ierr = PetscOptionsEnum("-ref_cell", "Use the reference cell", "ex40.c", DMPolytopeTypes, (PetscEnum) options->refCell, (PetscEnum *) &options->refCell, NULL);CHKERRQ(ierr);
30   ierr = PetscOptionsBool("-simplex", "Flag for simplices", "ex40.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
31   ierr = PetscOptionsEnd();
32   PetscFunctionReturn(0);
33 }
34 
35 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
36 {
37   PetscErrorCode ierr;
38 
39   PetscFunctionBegin;
40   if (ctx->refCell != DM_POLYTOPE_UNKNOWN) {
41     ierr = DMPlexCreateReferenceCellByType(comm, ctx->refCell, dm);CHKERRQ(ierr);
42   } else {
43     ierr = DMPlexCreateBoxMesh(comm, ctx->dim, ctx->simplex, NULL, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
44   }
45   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
46   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
47   PetscFunctionReturn(0);
48 }
49 
50 int main(int argc, char **argv)
51 {
52   DM             dm;
53   AppCtx         ctx;
54   PetscErrorCode ierr;
55 
56   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
57   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
58   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
59   ierr = DMDestroy(&dm);CHKERRQ(ierr);
60   ierr = PetscFinalize();
61   return ierr;
62 }
63 
64 /*TEST
65   test:
66     suffix: ref_tri
67     args: -ref_cell triangle -dm_refine 2 -dm_plex_check_all
68 
69   test:
70     suffix: box_tri
71     requires: triangle
72     nsize: {{1 3 5}}
73     args: -dm_plex_box_faces 3,3 -dm_refine 2 -dm_plex_check_all
74 
75   test:
76     suffix: ref_quad
77     args: -ref_cell quadrilateral -dm_refine 2 -dm_plex_check_all
78 
79   test:
80     suffix: box_quad
81     nsize: {{1 3 5}}
82     args: -dm_plex_box_faces 3,3 -simplex 0 -dm_refine 2 -dm_plex_check_all
83 
84   test:
85     suffix: ref_tet
86     args: -ref_cell tetrahedron -dm_refine 2 -dm_plex_check_all
87 
88   test:
89     suffix: box_tet
90     requires: ctetgen
91     nsize: {{1 3 5}}
92     args: -dim 3 -dm_plex_box_faces 3,3,3 -dm_refine 2 -dm_plex_check_all
93 
94   test:
95     suffix: ref_hex
96     args: -ref_cell hexahedron -simplex 0 -dm_refine 2 -dm_plex_check_all
97 
98   test:
99     suffix: box_hex
100     nsize: {{1 3 5}}
101     args: -dim 3 -dm_plex_box_faces 3,3,3 -simplex 0 -dm_refine 2 -dm_plex_check_all
102 
103   test:
104     suffix: ref_trip
105     args: -ref_cell triangular_prism -dm_refine 2 -dm_plex_check_all
106 
107   test:
108     suffix: ref_tquad
109     args: -ref_cell tensor_quad -dm_refine 2 -dm_plex_check_all
110 
111   test:
112     suffix: ref_ttrip
113     args: -ref_cell tensor_triangular_prism -dm_refine 2 -dm_plex_check_all
114 
115   test:
116     suffix: ref_tquadp
117     args: -ref_cell tensor_quadrilateral_prism -dm_refine 2 -dm_plex_check_all
118 
119   test:
120     suffix: ref_tri_tobox
121     args: -ref_cell triangle -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
122 
123   test:
124     suffix: box_tri_tobox
125     requires: triangle
126     nsize: {{1 3 5}}
127     args: -dm_plex_box_faces 3,3 -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
128 
129   test:
130     suffix: ref_tet_tobox
131     args: -ref_cell tetrahedron -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
132 
133   test:
134     suffix: box_tet_tobox
135     requires: ctetgen
136     nsize: {{1 3 5}}
137     args: -dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
138 
139   test:
140     suffix: ref_trip_tobox
141     args: -ref_cell triangular_prism -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
142 
143   test:
144     suffix: ref_ttrip_tobox
145     args: -ref_cell tensor_triangular_prism -dm_plex_cell_refiner tobox -dm_refine 2 -dm_plex_check_all
146 
147 TEST*/
148