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