xref: /petsc/src/dm/impls/plex/tests/ex66.c (revision d7d777e00e1db8154834d2ee4d21fc4a0c0a7f60)
1e7c9d70eSJeremy Theler static const char help[] = "Test for non-manifold interpolation";
2e7c9d70eSJeremy Theler 
3e7c9d70eSJeremy Theler #include <petscdmplex.h>
4e7c9d70eSJeremy Theler 
5e7c9d70eSJeremy Theler /*
6*b629e4b4SMatthew G. Knepley Test 0:
7*b629e4b4SMatthew G. Knepley 
8e7c9d70eSJeremy Theler        3-------------7
9e7c9d70eSJeremy Theler       /|            /|
10e7c9d70eSJeremy Theler      / |           / |
11e7c9d70eSJeremy Theler     /  |          /  |
12e7c9d70eSJeremy Theler    1-------------5   |
13e7c9d70eSJeremy Theler    |   |         |   |
14e7c9d70eSJeremy Theler    |   |         |   |
15e7c9d70eSJeremy Theler    |   |         |   |
16e7c9d70eSJeremy Theler    |   |         |   |
17e7c9d70eSJeremy Theler    z   4---------|---8
18e7c9d70eSJeremy Theler    ^  /          |  /
19e7c9d70eSJeremy Theler    | y           | /
20e7c9d70eSJeremy Theler    |/            |/
21e7c9d70eSJeremy Theler    2--->-x-------6-------------9
22*b629e4b4SMatthew G. Knepley 
23*b629e4b4SMatthew G. Knepley Test 1:
24*b629e4b4SMatthew G. Knepley 
25*b629e4b4SMatthew G. Knepley        3-------------7
26*b629e4b4SMatthew G. Knepley       /|            /|
27*b629e4b4SMatthew G. Knepley      / |           / |
28*b629e4b4SMatthew G. Knepley     /  |          /  |
29*b629e4b4SMatthew G. Knepley    1-------------5   |
30*b629e4b4SMatthew G. Knepley    |   |         |   |
31*b629e4b4SMatthew G. Knepley    |   |         |   |
32*b629e4b4SMatthew G. Knepley    |   |         |   |
33*b629e4b4SMatthew G. Knepley    |   |         |   |
34*b629e4b4SMatthew G. Knepley    z   4---------|---8
35*b629e4b4SMatthew G. Knepley    ^  /          |  / \
36*b629e4b4SMatthew G. Knepley    | y           | /   \
37*b629e4b4SMatthew G. Knepley    |/            |/     \
38*b629e4b4SMatthew G. Knepley    2--->-x-------6-------9
39e7c9d70eSJeremy Theler */
main(int argc,char ** argv)40e7c9d70eSJeremy Theler int main(int argc, char **argv)
41e7c9d70eSJeremy Theler {
42e7c9d70eSJeremy Theler   DM       dm, idm;
43e7c9d70eSJeremy Theler   DMLabel  ctLabel;
44*b629e4b4SMatthew G. Knepley   PetscInt testNum = 0;
45e7c9d70eSJeremy Theler 
46*b629e4b4SMatthew G. Knepley   PetscFunctionBeginUser;
47*b629e4b4SMatthew G. Knepley   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
48*b629e4b4SMatthew G. Knepley   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Non-Manifold Options", "ex66");
49*b629e4b4SMatthew G. Knepley   PetscCall(PetscOptionsInt("-test_num", "Test number", "", testNum, &testNum, NULL));
50*b629e4b4SMatthew G. Knepley   PetscOptionsEnd();
51*b629e4b4SMatthew G. Knepley 
52*b629e4b4SMatthew G. Knepley   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
53*b629e4b4SMatthew G. Knepley   PetscCall(PetscObjectSetName((PetscObject)dm, "cubeline-fromdag"));
54*b629e4b4SMatthew G. Knepley   PetscCall(DMSetType(dm, DMPLEX));
55*b629e4b4SMatthew G. Knepley   PetscCall(DMSetDimension(dm, 3));
56*b629e4b4SMatthew G. Knepley 
57*b629e4b4SMatthew G. Knepley   switch (testNum) {
58*b629e4b4SMatthew G. Knepley   case 0: {
59e7c9d70eSJeremy Theler     // 9 vertices
60e7c9d70eSJeremy Theler     // 1 edge
61e7c9d70eSJeremy Theler     // 0 faces
62e7c9d70eSJeremy Theler     // 1 volume
63e7c9d70eSJeremy Theler     PetscInt num_points[4] = {9, 1, 0, 1};
64e7c9d70eSJeremy Theler 
65e7c9d70eSJeremy Theler     // point 0 = hexahedron (defined by 8 vertices)
66e7c9d70eSJeremy Theler     // points 1-9 = vertices
67e7c9d70eSJeremy Theler     // point 10 = edged (defined by 2 vertices)
68e7c9d70eSJeremy Theler     PetscInt cone_size[11] = {8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 2};
69e7c9d70eSJeremy Theler 
70e7c9d70eSJeremy Theler     // hexahedron defined by points
71e7c9d70eSJeremy Theler     PetscInt    cones[11]             = {3, 4, 2, 1, 7, 5, 6, 8, 6, 9};
72e7c9d70eSJeremy Theler     PetscInt    cone_orientations[11] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
73e7c9d70eSJeremy Theler     PetscScalar vertex_coords[3 * 9]  = {0, 0, 1, 0, 0, 0, 0, 1, 1, 0, 1, 0, 1, 0, 1, 1, 0, 0, 1, 1, 1, 1, 1, 0, 2, 0, 0};
74e7c9d70eSJeremy Theler 
75*b629e4b4SMatthew G. Knepley     PetscCall(DMPlexCreateFromDAG(dm, 3, num_points, cone_size, cones, cone_orientations, vertex_coords));
76*b629e4b4SMatthew G. Knepley   } break;
77*b629e4b4SMatthew G. Knepley   case 1: {
78*b629e4b4SMatthew G. Knepley     // 9 vertices
79*b629e4b4SMatthew G. Knepley     // 0 edges
80*b629e4b4SMatthew G. Knepley     // 1 face
81*b629e4b4SMatthew G. Knepley     // 1 volume
82*b629e4b4SMatthew G. Knepley     PetscInt num_points[4] = {9, 0, 1, 1};
83e7c9d70eSJeremy Theler 
84*b629e4b4SMatthew G. Knepley     // point 0 = hexahedron (defined by 8 vertices)
85*b629e4b4SMatthew G. Knepley     // points 1-9 = vertices
86*b629e4b4SMatthew G. Knepley     // point 10 = triangle (defined by 3 vertices)
87*b629e4b4SMatthew G. Knepley     PetscInt cone_size[11] = {8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 3};
88*b629e4b4SMatthew G. Knepley 
89*b629e4b4SMatthew G. Knepley     // hexahedron defined by 8 points (point 0)
90*b629e4b4SMatthew G. Knepley     PetscInt    cones[11]             = {3, 4, 2, 1, 7, 5, 6, 8, 6, 9, 8};
91*b629e4b4SMatthew G. Knepley     PetscInt    cone_orientations[11] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
92*b629e4b4SMatthew G. Knepley     PetscScalar vertex_coords[3 * 9]  = {0, 0, 1,  // point 1
93*b629e4b4SMatthew G. Knepley                                          0, 0, 0,  // point 2
94*b629e4b4SMatthew G. Knepley                                          0, 1, 1,  // point 3
95*b629e4b4SMatthew G. Knepley                                          0, 1, 0,  // point 4
96*b629e4b4SMatthew G. Knepley                                          1, 0, 1,  // point 5
97*b629e4b4SMatthew G. Knepley                                          1, 0, 0,  // point 6
98*b629e4b4SMatthew G. Knepley                                          1, 1, 1,  // point 7
99*b629e4b4SMatthew G. Knepley                                          1, 1, 0,  // point 8
100*b629e4b4SMatthew G. Knepley                                          2, 0, 0}; // point 9
101e7c9d70eSJeremy Theler 
102e7c9d70eSJeremy Theler     PetscCall(DMPlexCreateFromDAG(dm, 3, num_points, cone_size, cones, cone_orientations, vertex_coords));
103*b629e4b4SMatthew G. Knepley   } break;
104*b629e4b4SMatthew G. Knepley   default:
105*b629e4b4SMatthew G. Knepley     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Invalid test number %" PetscInt_FMT, testNum);
106*b629e4b4SMatthew G. Knepley   }
107e7c9d70eSJeremy Theler   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
108e7c9d70eSJeremy Theler 
109e7c9d70eSJeremy Theler   // TODO: make it work with a DM made from a msh file
110e7c9d70eSJeremy Theler   // PetscCall(DMPlexCreateGmshFromFile(PETSC_COMM_WORLD, "cube-line.msh", PETSC_FALSE, &dm));
111e7c9d70eSJeremy Theler   // PetscCall(PetscObjectSetName((PetscObject)dm, "msh"));
112e7c9d70eSJeremy Theler 
113e7c9d70eSJeremy Theler   // Must set cell types
114e7c9d70eSJeremy Theler   PetscCall(DMPlexGetCellTypeLabel(dm, &ctLabel));
115*b629e4b4SMatthew G. Knepley   switch (testNum) {
116*b629e4b4SMatthew G. Knepley   case 0:
117e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 0, DM_POLYTOPE_HEXAHEDRON));
118e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 1, DM_POLYTOPE_POINT));
119e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 2, DM_POLYTOPE_POINT));
120e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 3, DM_POLYTOPE_POINT));
121e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 4, DM_POLYTOPE_POINT));
122e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 5, DM_POLYTOPE_POINT));
123e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 6, DM_POLYTOPE_POINT));
124e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 7, DM_POLYTOPE_POINT));
125e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 8, DM_POLYTOPE_POINT));
126e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 9, DM_POLYTOPE_POINT));
127e7c9d70eSJeremy Theler     PetscCall(DMLabelSetValue(ctLabel, 10, DM_POLYTOPE_SEGMENT));
128*b629e4b4SMatthew G. Knepley     break;
129*b629e4b4SMatthew G. Knepley   case 1:
130*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 0, DM_POLYTOPE_HEXAHEDRON));
131*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 1, DM_POLYTOPE_POINT));
132*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 2, DM_POLYTOPE_POINT));
133*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 3, DM_POLYTOPE_POINT));
134*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 4, DM_POLYTOPE_POINT));
135*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 5, DM_POLYTOPE_POINT));
136*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 6, DM_POLYTOPE_POINT));
137*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 7, DM_POLYTOPE_POINT));
138*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 8, DM_POLYTOPE_POINT));
139*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 9, DM_POLYTOPE_POINT));
140*b629e4b4SMatthew G. Knepley     PetscCall(DMLabelSetValue(ctLabel, 10, DM_POLYTOPE_TRIANGLE));
141*b629e4b4SMatthew G. Knepley     break;
142*b629e4b4SMatthew G. Knepley   default:
143*b629e4b4SMatthew G. Knepley     SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Invalid test number %" PetscInt_FMT, testNum);
144*b629e4b4SMatthew G. Knepley   }
145e7c9d70eSJeremy Theler 
146e7c9d70eSJeremy Theler   // interpolate (make sure to use -interp_dm_plex_stratify_celltype)
147e7c9d70eSJeremy Theler   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "interp_"));
148e7c9d70eSJeremy Theler   PetscCall(DMPlexInterpolate(dm, &idm));
149e7c9d70eSJeremy Theler   PetscCall(DMDestroy(&dm));
150e7c9d70eSJeremy Theler   dm = idm;
151e7c9d70eSJeremy Theler 
152e7c9d70eSJeremy Theler   PetscCall(DMSetFromOptions(dm));
153e7c9d70eSJeremy Theler   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
154e7c9d70eSJeremy Theler 
155e7c9d70eSJeremy Theler   PetscCall(DMDestroy(&dm));
156e7c9d70eSJeremy Theler   PetscCall(PetscFinalize());
157e7c9d70eSJeremy Theler   return 0;
158e7c9d70eSJeremy Theler }
159e7c9d70eSJeremy Theler 
160e7c9d70eSJeremy Theler /*TEST
161*b629e4b4SMatthew G. Knepley 
162*b629e4b4SMatthew G. Knepley   testset:
163*b629e4b4SMatthew G. Knepley     args: -interp_dm_plex_stratify_celltype -dm_view ::ascii_info_detail -interp_dm_view ::ascii_info_detail
164*b629e4b4SMatthew G. Knepley 
165e7c9d70eSJeremy Theler     test:
166e7c9d70eSJeremy Theler       suffix: 0
167*b629e4b4SMatthew G. Knepley     test:
168*b629e4b4SMatthew G. Knepley       suffix: 1
169*b629e4b4SMatthew G. Knepley       args: -test_num 1
170e7c9d70eSJeremy Theler 
171e7c9d70eSJeremy Theler TEST*/
172