xref: /petsc/src/dm/impls/plex/tests/ex4.c (revision 54fcfd0cd9e0b528eddf17b6a03728eb052cdf36)
1*54fcfd0cSMatthew G. Knepley static char help[] = "Tests for refinement of meshes created by hand\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   PetscInt  debug;          /* The debugging level */
7c4762a1bSJed Brown   PetscInt  dim;            /* The topological mesh dimension */
8c4762a1bSJed Brown   PetscBool cellHybrid;     /* Use a hybrid mesh */
9c4762a1bSJed Brown   PetscBool cellSimplex;    /* Use simplices or hexes */
10c4762a1bSJed Brown   PetscBool testPartition;  /* Use a fixed partitioning for testing */
11c4762a1bSJed Brown   PetscInt  testNum;        /* The particular mesh to test */
12c4762a1bSJed Brown   PetscBool uninterpolate;  /* Uninterpolate the mesh at the end */
13c4762a1bSJed Brown   PetscBool reinterpolate;  /* Reinterpolate the mesh at the end */
14c4762a1bSJed Brown } AppCtx;
15c4762a1bSJed Brown 
16c4762a1bSJed Brown PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   PetscErrorCode ierr;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   PetscFunctionBegin;
21c4762a1bSJed Brown   options->debug          = 0;
22c4762a1bSJed Brown   options->dim            = 2;
23c4762a1bSJed Brown   options->cellHybrid     = PETSC_TRUE;
24c4762a1bSJed Brown   options->cellSimplex    = PETSC_TRUE;
25c4762a1bSJed Brown   options->testPartition  = PETSC_TRUE;
26c4762a1bSJed Brown   options->testNum        = 0;
27c4762a1bSJed Brown   options->uninterpolate  = PETSC_FALSE;
28c4762a1bSJed Brown   options->reinterpolate  = PETSC_FALSE;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL,0);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
33c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL);CHKERRQ(ierr);
34c4762a1bSJed Brown   ierr = PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL);CHKERRQ(ierr);
36c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL,0);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = PetscOptionsEnd();
40c4762a1bSJed Brown   PetscFunctionReturn(0);
41c4762a1bSJed Brown }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown /* Two segments
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   2-------0-------3-------1-------4
46c4762a1bSJed Brown 
47c4762a1bSJed Brown become
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   4---0---7---1---5---2---8---3---6
50c4762a1bSJed Brown 
51c4762a1bSJed Brown */
52c4762a1bSJed Brown PetscErrorCode CreateSimplex_1D(MPI_Comm comm, DM *dm)
53c4762a1bSJed Brown {
54c4762a1bSJed Brown   PetscInt       depth = 1;
55c4762a1bSJed Brown   PetscMPIInt    rank;
56c4762a1bSJed Brown   PetscErrorCode ierr;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   PetscFunctionBegin;
59c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
60c4762a1bSJed Brown   if (!rank) {
61c4762a1bSJed Brown     PetscInt    numPoints[2]         = {3, 2};
62c4762a1bSJed Brown     PetscInt    coneSize[5]          = {2, 2, 0, 0, 0};
63c4762a1bSJed Brown     PetscInt    cones[4]             = {2, 3,  3, 4};
64c4762a1bSJed Brown     PetscInt    coneOrientations[16] = {0, 0,  0, 0};
65c4762a1bSJed Brown     PetscScalar vertexCoords[3]      = {-1.0, 0.0, 1.0};
66c4762a1bSJed Brown 
67c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
68c4762a1bSJed Brown   } else {
69c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
70c4762a1bSJed Brown 
71c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /* Two triangles
78c4762a1bSJed Brown         4
79c4762a1bSJed Brown       / | \
80c4762a1bSJed Brown      8  |  10
81c4762a1bSJed Brown     /   |   \
82c4762a1bSJed Brown    2  0 7  1 5
83c4762a1bSJed Brown     \   |   /
84c4762a1bSJed Brown      6  |  9
85c4762a1bSJed Brown       \ | /
86c4762a1bSJed Brown         3
87c4762a1bSJed Brown 
88c4762a1bSJed Brown Becomes
89c4762a1bSJed Brown            10
90c4762a1bSJed Brown           / | \
91c4762a1bSJed Brown         21  |  26
92c4762a1bSJed Brown         /   |   \
93c4762a1bSJed Brown       14 2 20 4  16
94c4762a1bSJed Brown       /|\   |   /|\
95c4762a1bSJed Brown     22 | 28 | 32 | 25
96c4762a1bSJed Brown     /  |  \ | /  | 6\
97c4762a1bSJed Brown    8  29 3 13  7 31  11
98c4762a1bSJed Brown     \0 |  / | \  |  /
99c4762a1bSJed Brown     17 | 27 | 30 | 24
100c4762a1bSJed Brown       \|/   |   \|/
101c4762a1bSJed Brown       12 1 19 5  15
102c4762a1bSJed Brown         \   |   /
103c4762a1bSJed Brown         18  |  23
104c4762a1bSJed Brown           \ | /
105c4762a1bSJed Brown             9
106c4762a1bSJed Brown */
107c4762a1bSJed Brown PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm)
108c4762a1bSJed Brown {
109c4762a1bSJed Brown   PetscInt       depth = 2;
110c4762a1bSJed Brown   PetscMPIInt    rank;
111c4762a1bSJed Brown   PetscErrorCode ierr;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBegin;
114c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
115c4762a1bSJed Brown   if (!rank) {
116c4762a1bSJed Brown     PetscInt    numPoints[3]         = {4, 5, 2};
117c4762a1bSJed Brown     PetscInt    coneSize[11]         = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2};
118c4762a1bSJed Brown     PetscInt    cones[16]            = {6, 7, 8,  7, 9, 10,  2, 3,  3, 4,  4, 2,  3, 5,  5, 4};
119c4762a1bSJed Brown     PetscInt    coneOrientations[16] = {0, 0, 0, -2, 0,  0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0};
120c4762a1bSJed Brown     PetscScalar vertexCoords[8]      = {-0.5, 0.0,  0.0, -0.5,  0.0, 0.5,  0.5, 0.0};
121c4762a1bSJed Brown 
122c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
123c4762a1bSJed Brown   } else {
124c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
127c4762a1bSJed Brown   }
128c4762a1bSJed Brown   PetscFunctionReturn(0);
129c4762a1bSJed Brown }
130c4762a1bSJed Brown 
131c4762a1bSJed Brown /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges
132c4762a1bSJed Brown         5--16--8
133c4762a1bSJed Brown       / |      | \
134c4762a1bSJed Brown     11  |      |  12
135c4762a1bSJed Brown     /   |      |   \
136c4762a1bSJed Brown    3  0 10  2 14 1  6
137c4762a1bSJed Brown     \   |      |   /
138c4762a1bSJed Brown      9  |      |  13
139c4762a1bSJed Brown       \ |      | /
140c4762a1bSJed Brown         4--15--7
141c4762a1bSJed Brown */
142c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
143c4762a1bSJed Brown {
144c4762a1bSJed Brown   DM             idm, hdm = NULL;
145c4762a1bSJed Brown   DMLabel        faultLabel, hybridLabel;
146c4762a1bSJed Brown   PetscInt       p;
147c4762a1bSJed Brown   PetscMPIInt    rank;
148c4762a1bSJed Brown   PetscErrorCode ierr;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBegin;
151c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
152c4762a1bSJed Brown   if (!rank) {
153c4762a1bSJed Brown     switch (testNum) {
154c4762a1bSJed Brown     case 0:
155c4762a1bSJed Brown     {
156c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
157c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
158c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4,  5, 4, 3};
159c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0,  0, 0, 0};
160c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  1.0, 0.5};
161c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
162c4762a1bSJed Brown 
163c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
164c4762a1bSJed Brown       for (p = 0; p < 2; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
165c4762a1bSJed Brown     }
166c4762a1bSJed Brown     break;
167c4762a1bSJed Brown     case 1:
168c4762a1bSJed Brown     {
169c4762a1bSJed Brown       PetscInt    numPoints[2]         = {5, 4};
170c4762a1bSJed Brown       PetscInt    coneSize[9]          = {3, 3, 3, 3, 0, 0, 0, 0, 0};
171c4762a1bSJed Brown       PetscInt    cones[12]            = {4, 5, 6,  6, 7, 4,  6, 5, 8,  6, 8, 7};
172c4762a1bSJed Brown       PetscInt    coneOrientations[12] = {0, 0, 0,  0, 0, 0,  0, 0, 0,  0, 0, 0};
173c4762a1bSJed Brown       PetscScalar vertexCoords[10]     = {-1.0, 0.0,  0.0, -1.0,  0.0, 0.0,  0.0, 1.0,  1.0, 0.0};
174c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {5, 6, 7};
175c4762a1bSJed Brown 
176c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
177c4762a1bSJed Brown       for (p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown     break;
180c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
181c4762a1bSJed Brown     }
182c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
185c4762a1bSJed Brown     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
186c4762a1bSJed Brown     ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr);
187c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr);
188c4762a1bSJed Brown     ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr);
189c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
190c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
191c4762a1bSJed Brown     *dm  = hdm;
192c4762a1bSJed Brown   } else {
193c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
196c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
197*54fcfd0cSMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
198*54fcfd0cSMatthew G. Knepley     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
199*54fcfd0cSMatthew G. Knepley     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
200c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr);
201c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
202c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
203c4762a1bSJed Brown     *dm  = hdm;
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown /* Two quadrilaterals
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   5----10-----4----14-----7
211c4762a1bSJed Brown   |           |           |
212c4762a1bSJed Brown   |           |           |
213c4762a1bSJed Brown   |           |           |
214c4762a1bSJed Brown  11     0     9     1     13
215c4762a1bSJed Brown   |           |           |
216c4762a1bSJed Brown   |           |           |
217c4762a1bSJed Brown   |           |           |
218c4762a1bSJed Brown   2-----8-----3----12-----6
219c4762a1bSJed Brown */
220c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   PetscInt       depth = 2;
223c4762a1bSJed Brown   PetscMPIInt    rank;
224c4762a1bSJed Brown   PetscErrorCode ierr;
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   PetscFunctionBegin;
227c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
228c4762a1bSJed Brown   if (!rank) {
229c4762a1bSJed Brown     PetscInt    numPoints[3]         = {6, 7, 2};
230c4762a1bSJed Brown     PetscInt    coneSize[15]         = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2};
231c4762a1bSJed Brown     PetscInt    cones[22]            = {8, 9, 10, 11,  12, 13, 14,  9,  2, 3,  3, 4,  4, 5,  5, 2,  3, 6,  6, 7,  7, 4};
232c4762a1bSJed Brown     PetscInt    coneOrientations[22] = {0, 0,  0,  0,   0,  0,  0, -2,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0};
233c4762a1bSJed Brown     PetscScalar vertexCoords[12]     = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
234c4762a1bSJed Brown 
235c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
236c4762a1bSJed Brown   } else {
237c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
238c4762a1bSJed Brown 
239c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
240c4762a1bSJed Brown   }
241c4762a1bSJed Brown   PetscFunctionReturn(0);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
245c4762a1bSJed Brown {
246c4762a1bSJed Brown   DM             idm, hdm = NULL;
247c4762a1bSJed Brown   DMLabel        faultLabel, hybridLabel;
248c4762a1bSJed Brown   PetscInt       p;
249c4762a1bSJed Brown   PetscMPIInt    rank;
250c4762a1bSJed Brown   PetscErrorCode ierr;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBegin;
253c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
254c4762a1bSJed Brown   if (!rank) {
255c4762a1bSJed Brown     PetscInt    numPoints[2]        = {6, 2};
256c4762a1bSJed Brown     PetscInt    coneSize[8]         = {4, 4, 0, 0, 0, 0, 0, 0};
257c4762a1bSJed Brown     PetscInt    cones[8]            = {2, 3, 4, 5,  3, 6, 7, 4,};
258c4762a1bSJed Brown     PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
259c4762a1bSJed Brown     PetscScalar vertexCoords[12]    = {-1.0, -0.5,  0.0, -0.5,  0.0, 0.5,  -1.0, 0.5,  1.0, -0.5,  1.0, 0.5};
260c4762a1bSJed Brown     PetscInt    faultPoints[2]      = {3, 4};
261c4762a1bSJed Brown 
262c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
263c4762a1bSJed Brown     for(p = 0; p < 2; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
264c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
265c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
266c4762a1bSJed Brown     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
267c4762a1bSJed Brown     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
268c4762a1bSJed Brown     ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr);
269c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr);
270c4762a1bSJed Brown     ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr);
271c4762a1bSJed Brown   } else {
272c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
273c4762a1bSJed Brown 
274c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
275c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
276*54fcfd0cSMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
277*54fcfd0cSMatthew G. Knepley     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
278*54fcfd0cSMatthew G. Knepley     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
279c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr);
280c4762a1bSJed Brown   }
281c4762a1bSJed Brown   ierr = DMDestroy(&idm);CHKERRQ(ierr);
282c4762a1bSJed Brown   ierr = DMDestroy(dm);CHKERRQ(ierr);
283c4762a1bSJed Brown   *dm  = hdm;
284c4762a1bSJed Brown   PetscFunctionReturn(0);
285c4762a1bSJed Brown }
286c4762a1bSJed Brown 
287c4762a1bSJed Brown /* Two tetrahedrons
288c4762a1bSJed Brown 
289c4762a1bSJed Brown  cell   5          5______    cell
290c4762a1bSJed Brown  0    / | \        |\      \     1
291c4762a1bSJed Brown     17  |  18      | 18 13  21
292c4762a1bSJed Brown     /8 19 10\     19  \      \
293c4762a1bSJed Brown    2-14-|----4     |   4--22--6
294c4762a1bSJed Brown     \ 9 | 7 /      |10 /      /
295c4762a1bSJed Brown     16  |  15      | 15  12 20
296c4762a1bSJed Brown       \ | /        |/      /
297c4762a1bSJed Brown         3          3------
298c4762a1bSJed Brown */
299c4762a1bSJed Brown PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
300c4762a1bSJed Brown {
301c4762a1bSJed Brown   DM             idm;
302c4762a1bSJed Brown   PetscInt       depth = 3;
303c4762a1bSJed Brown   PetscMPIInt    rank;
304c4762a1bSJed Brown   PetscErrorCode ierr;
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   PetscFunctionBegin;
307c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
308c4762a1bSJed Brown   if (!rank) {
309c4762a1bSJed Brown     switch (testNum) {
310c4762a1bSJed Brown     case 0:
311c4762a1bSJed Brown     {
312c4762a1bSJed Brown       PetscInt    numPoints[4]         = {5, 9, 7, 2};
313c4762a1bSJed Brown       PetscInt    coneSize[23]         = {4, 4, 0, 0, 0, 0, 0, 3, 3, 3, 3, 3, 3, 3, 2, 2, 2, 2, 2, 2, 2, 2, 2};
314c4762a1bSJed Brown       PetscInt    cones[47]            = { 7,  8,  9, 10,  10, 11, 12, 13,  14, 15, 16,  17, 18, 14,  16, 19, 17,  15, 18, 19,  20, 21, 19,  15, 22, 20,  18, 21, 22,  2, 4,  4, 3,  3, 2,  2, 5,  5, 4,  3, 5,  3, 6,  6, 5,  4, 6};
315c4762a1bSJed Brown       PetscInt    coneOrientations[47] = { 0,  0,  0,  0,  -3,  0,  0,  0,   0,  0,  0,   0,  0, -2,  -2,  0, -2,  -2, -2, -2,   0,  0, -2,  -2,  0, -2,  -2, -2, -2,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0,  0, 0};
316c4762a1bSJed Brown       PetscScalar vertexCoords[15]     = {0.0, 0.0, -0.5,  0.0, -0.5, 0.0,  1.0, 0.0, 0.0,  0.0, 0.5, 0.0,  0.0, 0.0, 0.5};
317c4762a1bSJed Brown 
318c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
319c4762a1bSJed Brown     }
320c4762a1bSJed Brown     break;
321c4762a1bSJed Brown     case 1:
322c4762a1bSJed Brown     {
323c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
324c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
325c4762a1bSJed Brown       PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
326c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
327c4762a1bSJed Brown       PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
328c4762a1bSJed Brown 
329c4762a1bSJed Brown       depth = 1;
330c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
331c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
332c4762a1bSJed Brown       ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr);
333c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
334c4762a1bSJed Brown       *dm  = idm;
335c4762a1bSJed Brown     }
336c4762a1bSJed Brown     break;
337c4762a1bSJed Brown     case 2:
338c4762a1bSJed Brown     {
339c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 1};
340c4762a1bSJed Brown       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
341c4762a1bSJed Brown       PetscInt    cones[4]            = {2, 3, 4, 1};
342c4762a1bSJed Brown       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
343c4762a1bSJed Brown       PetscScalar vertexCoords[12]    = {0.0, 0.0, 0.0,  1.0, 0.0, 0.0,  0.0, 1.0, 0.0,  0.0, 0.0, 1.0};
344c4762a1bSJed Brown 
345c4762a1bSJed Brown       depth = 1;
346c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
347c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
348c4762a1bSJed Brown       ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr);
349c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
350c4762a1bSJed Brown       *dm  = idm;
351c4762a1bSJed Brown     }
352c4762a1bSJed Brown     break;
353c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
354c4762a1bSJed Brown     }
355c4762a1bSJed Brown   } else {
356c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
357c4762a1bSJed Brown 
358c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
359c4762a1bSJed Brown     switch (testNum) {
360c4762a1bSJed Brown     case 1:
361c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
362c4762a1bSJed Brown       ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr);
363c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
364c4762a1bSJed Brown       *dm  = idm;
365c4762a1bSJed Brown       break;
366c4762a1bSJed Brown     }
367c4762a1bSJed Brown   }
368c4762a1bSJed Brown   PetscFunctionReturn(0);
369c4762a1bSJed Brown }
370c4762a1bSJed Brown 
371c4762a1bSJed Brown /* Two tetrahedrons separated by a zero-volume cell with 6 vertices
372c4762a1bSJed Brown 
373c4762a1bSJed Brown  cell   6 ___33___10______    cell
374c4762a1bSJed Brown  0    / | \        |\      \     1
375c4762a1bSJed Brown     21  |  23      | 29     27
376c4762a1bSJed Brown     /12 24 14\    30  \      \
377c4762a1bSJed Brown    3-20-|----5--32-|---9--26--7
378c4762a1bSJed Brown     \ 13| 11/      |18 /      /
379c4762a1bSJed Brown     19  |  22      | 28     25
380c4762a1bSJed Brown       \ | /        |/      /
381c4762a1bSJed Brown         4----31----8------
382c4762a1bSJed Brown          cell 2
383c4762a1bSJed Brown */
384c4762a1bSJed Brown PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
385c4762a1bSJed Brown {
386c4762a1bSJed Brown   DM             idm, hdm = NULL;
387c4762a1bSJed Brown   DMLabel        faultLabel, hybridLabel;
388c4762a1bSJed Brown   PetscInt       p;
389c4762a1bSJed Brown   PetscMPIInt    rank;
390c4762a1bSJed Brown   PetscErrorCode ierr;
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   PetscFunctionBegin;
393c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
394c4762a1bSJed Brown   if (!rank) {
395c4762a1bSJed Brown     switch (testNum) {
396c4762a1bSJed Brown     case 0:
397c4762a1bSJed Brown     {
398c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
399c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
400c4762a1bSJed Brown       PetscInt    cones[8]            = {4, 3, 5, 2,  5, 3, 4, 6};
401c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0,  0, 0, 0, 0};
402c4762a1bSJed Brown       PetscScalar vertexCoords[15]    = {-1.0, 0.0, 0.0,  0.0, -1.0, 0.0,  0.0, 0.0, 1.0,  0.0, 1.0, 0.0,  1.0, 0.0, 0.0};
403c4762a1bSJed Brown       PetscInt    faultPoints[3]      = {3, 4, 5};
404c4762a1bSJed Brown 
405c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
406c4762a1bSJed Brown       for(p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
407c4762a1bSJed Brown     }
408c4762a1bSJed Brown     break;
409c4762a1bSJed Brown     case 1:
410c4762a1bSJed Brown     {
411c4762a1bSJed Brown       /* Tets 0,3,5 and 1,2,4 */
412c4762a1bSJed Brown       PetscInt    numPoints[2]         = {9, 6};
413c4762a1bSJed Brown       PetscInt    coneSize[15]         = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0};
414*54fcfd0cSMatthew G. Knepley       PetscInt    cones[24]            = { 7,  8,  9, 6,  11, 13,  9, 14,  10, 11, 13, 9,
415*54fcfd0cSMatthew G. Knepley                                            9, 10, 11, 7,   9, 14, 13, 12,   7,  8, 11, 9};
416c4762a1bSJed Brown       PetscInt    coneOrientations[24] = { 0, 0,  0, 0,   0,  0,  0,  0,   0,  0,  0, 0,
417c4762a1bSJed Brown                                            0, 0,  0, 0,   0,  0,  0,  0,   0,  0,  0, 0};
418c4762a1bSJed Brown       PetscScalar vertexCoords[27]     = {-2.0, -1.0,  0.0,  -2.0,  0.0,  0.0,  -2.0,  0.0,  1.0,
419c4762a1bSJed Brown                                            0.0, -1.0,  0.0,   0.0,  0.0,  0.0,   0.0,  0.0,  1.0,
420c4762a1bSJed Brown                                            2.0, -1.0,  0.0,   2.0,  0.0,  0.0,   2.0,  0.0,  1.0};
421c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {9, 10, 11};
422c4762a1bSJed Brown 
423c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
424c4762a1bSJed Brown       for(p = 0; p < 3; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
425c4762a1bSJed Brown     }
426c4762a1bSJed Brown     break;
427c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
428c4762a1bSJed Brown     }
429c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
430c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
431c4762a1bSJed Brown     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
432c4762a1bSJed Brown     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
433c4762a1bSJed Brown     ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr);
434c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, &hybridLabel, NULL, NULL, &hdm);CHKERRQ(ierr);
435c4762a1bSJed Brown     ierr = DMLabelDestroy(&hybridLabel);CHKERRQ(ierr);
436c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
437c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
438c4762a1bSJed Brown     *dm  = hdm;
439c4762a1bSJed Brown   } else {
440c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
441c4762a1bSJed Brown 
442c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
443c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
444*54fcfd0cSMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
445*54fcfd0cSMatthew G. Knepley     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
446*54fcfd0cSMatthew G. Knepley     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
447c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr);
448c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
449c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
450c4762a1bSJed Brown     *dm  = hdm;
451c4762a1bSJed Brown   }
452c4762a1bSJed Brown   PetscFunctionReturn(0);
453c4762a1bSJed Brown }
454c4762a1bSJed Brown 
455c4762a1bSJed Brown PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
456c4762a1bSJed Brown {
457c4762a1bSJed Brown   DM             idm;
458c4762a1bSJed Brown   PetscMPIInt    rank;
459c4762a1bSJed Brown   PetscErrorCode ierr;
460c4762a1bSJed Brown 
461c4762a1bSJed Brown   PetscFunctionBegin;
462c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
463c4762a1bSJed Brown   if (!rank) {
464c4762a1bSJed Brown     switch (testNum) {
465c4762a1bSJed Brown     case 0:
466c4762a1bSJed Brown     {
467c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
468c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
469c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
470c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
471c4762a1bSJed Brown       PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
472c4762a1bSJed Brown                                           -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
473c4762a1bSJed Brown                                           1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
474c4762a1bSJed Brown 
475c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
476c4762a1bSJed Brown     }
477c4762a1bSJed Brown     break;
478c4762a1bSJed Brown     case 1:
479c4762a1bSJed Brown     {
480c4762a1bSJed Brown       PetscInt    numPoints[2]        = {8, 1};
481c4762a1bSJed Brown       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
482c4762a1bSJed Brown       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
483c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
484c4762a1bSJed Brown       PetscScalar vertexCoords[24]    = {-1.0, -1.0, -1.0,  -1.0,  1.0, -1.0,  1.0,  1.0, -1.0,   1.0, -1.0, -1.0,
485c4762a1bSJed Brown                                          -1.0, -1.0,  1.0,   1.0, -1.0,  1.0,  1.0,  1.0,  1.0,  -1.0,  1.0,  1.0};
486c4762a1bSJed Brown 
487c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
488c4762a1bSJed Brown     }
489c4762a1bSJed Brown     break;
490c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
491c4762a1bSJed Brown     }
492c4762a1bSJed Brown   } else {
493c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
494c4762a1bSJed Brown 
495c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
496c4762a1bSJed Brown   }
497c4762a1bSJed Brown   ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
498c4762a1bSJed Brown   ierr = DMViewFromOptions(idm, NULL, "-in_dm_view");CHKERRQ(ierr);
499c4762a1bSJed Brown   ierr = DMDestroy(dm);CHKERRQ(ierr);
500c4762a1bSJed Brown   *dm  = idm;
501c4762a1bSJed Brown   PetscFunctionReturn(0);
502c4762a1bSJed Brown }
503c4762a1bSJed Brown 
504c4762a1bSJed Brown PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
505c4762a1bSJed Brown {
506c4762a1bSJed Brown   DM             idm, hdm = NULL;
507c4762a1bSJed Brown   DMLabel        faultLabel;
508c4762a1bSJed Brown   PetscInt       p;
509c4762a1bSJed Brown   PetscMPIInt    rank;
510c4762a1bSJed Brown   PetscErrorCode ierr;
511c4762a1bSJed Brown 
512c4762a1bSJed Brown   PetscFunctionBegin;
513c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
514c4762a1bSJed Brown   if (!rank) {
515c4762a1bSJed Brown     switch (testNum) {
516c4762a1bSJed Brown     case 0:
517c4762a1bSJed Brown     {
518c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
519c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
520c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9,  5, 4, 10, 11, 7, 12, 13, 8};
521c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0,  0,  0, 0,  0,  0, 0};
522c4762a1bSJed Brown       PetscScalar vertexCoords[36]     = {-1.0, -0.5, -0.5,  -1.0,  0.5, -0.5,  0.0,  0.5, -0.5,   0.0, -0.5, -0.5,
523c4762a1bSJed Brown                                           -1.0, -0.5,  0.5,   0.0, -0.5,  0.5,  0.0,  0.5,  0.5,  -1.0,  0.5,  0.5,
524c4762a1bSJed Brown                                           1.0,  0.5, -0.5,   1.0, -0.5, -0.5,  1.0, -0.5,  0.5,   1.0,  0.5,  0.5};
525c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {2, 3, 5, 6};
526c4762a1bSJed Brown 
527c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
528c4762a1bSJed Brown       for(p = 0; p < 4; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
529c4762a1bSJed Brown     }
530c4762a1bSJed Brown     break;
531c4762a1bSJed Brown     case 1:
532c4762a1bSJed Brown     {
533c4762a1bSJed Brown       PetscInt    numPoints[2]         = {30, 7};
534c4762a1bSJed Brown       PetscInt    coneSize[37]         = {8,8,8,8,8,8,8, 0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
535c4762a1bSJed Brown       PetscInt    cones[56]            = { 8, 21, 20,  7, 13, 12, 23, 24,
536c4762a1bSJed Brown                                           14, 15, 10,  9, 13,  8, 21, 24,
537c4762a1bSJed Brown                                           15, 16, 11, 10, 24, 21, 22, 25,
538c4762a1bSJed Brown                                           30, 29, 28, 21, 35, 24, 33, 34,
539c4762a1bSJed Brown                                           24, 21, 30, 35, 25, 36, 31, 22,
540c4762a1bSJed Brown                                           27, 20, 21, 28, 32, 33, 24, 23,
541c4762a1bSJed Brown                                           15, 24, 13, 14, 19, 18, 17, 26};
542c4762a1bSJed Brown       PetscInt    coneOrientations[56] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
543c4762a1bSJed Brown       PetscScalar vertexCoords[90]     = {-2.0, -2.0, -2.0,  -2.0, -1.0, -2.0,  -3.0,  0.0, -2.0,  -2.0,  1.0, -2.0,  -2.0,  2.0, -2.0,  -2.0, -2.0,  0.0,
544c4762a1bSJed Brown                                           -2.0, -1.0,  0.0,  -3.0,  0.0,  0.0,  -2.0,  1.0,  0.0,  -2.0,  2.0,  0.0,  -2.0, -1.0,  2.0,  -3.0,  0.0,  2.0,
545c4762a1bSJed Brown                                           -2.0,  1.0,  2.0,   0.0, -2.0, -2.0,   0.0,  0.0, -2.0,   0.0,  2.0, -2.0,   0.0, -2.0,  0.0,   0.0,  0.0,  0.0,
546c4762a1bSJed Brown                                            0.0,  2.0,  0.0,   0.0,  0.0,  2.0,   2.0, -2.0, -2.0,   2.0, -1.0, -2.0,   3.0,  0.0, -2.0,   2.0,  1.0, -2.0,
547c4762a1bSJed Brown                                            2.0,  2.0, -2.0,   2.0, -2.0,  0.0,   2.0, -1.0,  0.0,   3.0,  0.0,  0.0,   2.0,  1.0,  0.0,   2.0,  2.0,  0.0};
548c4762a1bSJed Brown       PetscInt    faultPoints[6]       = {20, 21, 22, 23, 24, 25};
549c4762a1bSJed Brown 
550c4762a1bSJed Brown       ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
551c4762a1bSJed Brown       for(p = 0; p < 6; ++p) {ierr = DMSetLabelValue(*dm, "fault", faultPoints[p], 1);CHKERRQ(ierr);}
552c4762a1bSJed Brown     }
553c4762a1bSJed Brown     break;
554c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %d", testNum);
555c4762a1bSJed Brown     }
556c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
557c4762a1bSJed Brown     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
558c4762a1bSJed Brown     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
559c4762a1bSJed Brown     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
560c4762a1bSJed Brown     ierr = DMGetLabel(*dm, "fault", &faultLabel);CHKERRQ(ierr);
561c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, faultLabel, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr);
562c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
563c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
564c4762a1bSJed Brown     *dm  = hdm;
565c4762a1bSJed Brown   } else {
566c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
567c4762a1bSJed Brown 
568c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL);CHKERRQ(ierr);
569c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
570*54fcfd0cSMatthew G. Knepley     ierr = PetscObjectSetOptionsPrefix((PetscObject) idm, "in_");CHKERRQ(ierr);
571*54fcfd0cSMatthew G. Knepley     ierr = DMSetFromOptions(idm);CHKERRQ(ierr);
572*54fcfd0cSMatthew G. Knepley     ierr = DMViewFromOptions(idm, NULL, "-dm_view");CHKERRQ(ierr);
573c4762a1bSJed Brown     ierr = DMPlexCreateHybridMesh(idm, NULL, NULL, NULL, NULL, NULL, &hdm);CHKERRQ(ierr);
574c4762a1bSJed Brown     ierr = DMDestroy(&idm);CHKERRQ(ierr);
575c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
576c4762a1bSJed Brown     *dm  = hdm;
577c4762a1bSJed Brown   }
578c4762a1bSJed Brown   PetscFunctionReturn(0);
579c4762a1bSJed Brown }
580c4762a1bSJed Brown 
581c4762a1bSJed Brown PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
582c4762a1bSJed Brown {
583c4762a1bSJed Brown   PetscInt       dim         = user->dim;
584c4762a1bSJed Brown   PetscBool      cellHybrid  = user->cellHybrid;
585c4762a1bSJed Brown   PetscBool      cellSimplex = user->cellSimplex;
586c4762a1bSJed Brown   PetscMPIInt    rank, size;
587c4762a1bSJed Brown   PetscErrorCode ierr;
588c4762a1bSJed Brown 
589c4762a1bSJed Brown   PetscFunctionBegin;
590c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
591c4762a1bSJed Brown   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
592c4762a1bSJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
593c4762a1bSJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
594c4762a1bSJed Brown   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
595c4762a1bSJed Brown   switch (dim) {
596c4762a1bSJed Brown   case 1:
597c4762a1bSJed Brown     if (cellHybrid) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %d", dim);
598c4762a1bSJed Brown     ierr = CreateSimplex_1D(comm, dm);CHKERRQ(ierr);
599c4762a1bSJed Brown     break;
600c4762a1bSJed Brown   case 2:
601c4762a1bSJed Brown     if (cellSimplex) {
602c4762a1bSJed Brown       if (cellHybrid) {
603c4762a1bSJed Brown         ierr = CreateSimplexHybrid_2D(comm, user->testNum, dm);CHKERRQ(ierr);
604c4762a1bSJed Brown       } else {
605c4762a1bSJed Brown         ierr = CreateSimplex_2D(comm, dm);CHKERRQ(ierr);
606c4762a1bSJed Brown       }
607c4762a1bSJed Brown     } else {
608c4762a1bSJed Brown       if (cellHybrid) {
609c4762a1bSJed Brown         ierr = CreateTensorProductHybrid_2D(comm, user->testNum, dm);CHKERRQ(ierr);
610c4762a1bSJed Brown       } else {
611c4762a1bSJed Brown         ierr = CreateTensorProduct_2D(comm, user->testNum, dm);CHKERRQ(ierr);
612c4762a1bSJed Brown       }
613c4762a1bSJed Brown     }
614c4762a1bSJed Brown     break;
615c4762a1bSJed Brown   case 3:
616c4762a1bSJed Brown     if (cellSimplex) {
617c4762a1bSJed Brown       if (cellHybrid) {
618c4762a1bSJed Brown         ierr = CreateSimplexHybrid_3D(comm, user->testNum, dm);CHKERRQ(ierr);
619c4762a1bSJed Brown       } else {
620c4762a1bSJed Brown         ierr = CreateSimplex_3D(comm, user->testNum, dm);CHKERRQ(ierr);
621c4762a1bSJed Brown       }
622c4762a1bSJed Brown     } else {
623c4762a1bSJed Brown       if (cellHybrid) {
624c4762a1bSJed Brown         ierr = CreateTensorProductHybrid_3D(comm, user->testNum, dm);CHKERRQ(ierr);
625c4762a1bSJed Brown       } else {
626c4762a1bSJed Brown         ierr = CreateTensorProduct_3D(comm, user->testNum, dm);CHKERRQ(ierr);
627c4762a1bSJed Brown       }
628c4762a1bSJed Brown     }
629c4762a1bSJed Brown     break;
630c4762a1bSJed Brown   default:
631c4762a1bSJed Brown     SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %d", dim);
632c4762a1bSJed Brown   }
633c4762a1bSJed Brown   if (user->testPartition && size > 1) {
634c4762a1bSJed Brown     PetscPartitioner part;
635c4762a1bSJed Brown     PetscInt  *sizes  = NULL;
636c4762a1bSJed Brown     PetscInt  *points = NULL;
637c4762a1bSJed Brown 
638c4762a1bSJed Brown     if (!rank) {
639c4762a1bSJed Brown       if (dim == 2 && cellSimplex && !cellHybrid && size == 2) {
640c4762a1bSJed Brown         switch (user->testNum) {
641c4762a1bSJed Brown         case 0: {
642c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 1};
643c4762a1bSJed Brown           PetscInt triPoints_p2[2] = {0, 1};
644c4762a1bSJed Brown 
645c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
646c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
647c4762a1bSJed Brown           ierr = PetscArraycpy(points, triPoints_p2, 2);CHKERRQ(ierr);break;}
648c4762a1bSJed Brown         default:
649c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular mesh on 2 procs", user->testNum);
650c4762a1bSJed Brown         }
651c4762a1bSJed Brown       } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) {
652c4762a1bSJed Brown         switch (user->testNum) {
653c4762a1bSJed Brown         case 0: {
654c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 2};
655c4762a1bSJed Brown           PetscInt triPoints_p2[3] = {0, 1, 2};
656c4762a1bSJed Brown 
657c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr);
658c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  triSizes_p2, 2);CHKERRQ(ierr);
659c4762a1bSJed Brown           ierr = PetscArraycpy(points, triPoints_p2, 3);CHKERRQ(ierr);break;}
660c4762a1bSJed Brown         default:
661c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for triangular hybrid mesh on 2 procs", user->testNum);
662c4762a1bSJed Brown         }
663c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) {
664c4762a1bSJed Brown         switch (user->testNum) {
665c4762a1bSJed Brown         case 0: {
666c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 1};
667c4762a1bSJed Brown           PetscInt quadPoints_p2[2] = {0, 1};
668c4762a1bSJed Brown 
669c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
670c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  quadSizes_p2, 2);CHKERRQ(ierr);
671c4762a1bSJed Brown           ierr = PetscArraycpy(points, quadPoints_p2, 2);CHKERRQ(ierr);break;}
672c4762a1bSJed Brown         default:
673c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral mesh on 2 procs", user->testNum);
674c4762a1bSJed Brown         }
675c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) {
676c4762a1bSJed Brown         switch (user->testNum) {
677c4762a1bSJed Brown         case 0: {
678c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 2};
679c4762a1bSJed Brown           PetscInt quadPoints_p2[3] = {0, 1, 2};
680c4762a1bSJed Brown 
681c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr);
682c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  quadSizes_p2, 2);CHKERRQ(ierr);
683c4762a1bSJed Brown           ierr = PetscArraycpy(points, quadPoints_p2, 3);CHKERRQ(ierr);break;}
684c4762a1bSJed Brown         default:
685c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for quadrilateral hybrid mesh on 2 procs", user->testNum);
686c4762a1bSJed Brown         }
687c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) {
688c4762a1bSJed Brown         switch (user->testNum) {
689c4762a1bSJed Brown         case 0: {
690c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 1};
691c4762a1bSJed Brown           PetscInt tetPoints_p2[2] = {0, 1};
692c4762a1bSJed Brown 
693c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
694c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  tetSizes_p2, 2);CHKERRQ(ierr);
695c4762a1bSJed Brown           ierr = PetscArraycpy(points, tetPoints_p2, 2);CHKERRQ(ierr);break;}
696c4762a1bSJed Brown         case 1: {
697c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 1};
698c4762a1bSJed Brown           PetscInt tetPoints_p2[2] = {0, 1};
699c4762a1bSJed Brown 
700c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
701c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  tetSizes_p2, 2);CHKERRQ(ierr);
702c4762a1bSJed Brown           ierr = PetscArraycpy(points, tetPoints_p2, 2);CHKERRQ(ierr);break;}
703c4762a1bSJed Brown         default:
704c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral mesh on 2 procs", user->testNum);
705c4762a1bSJed Brown         }
706c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) {
707c4762a1bSJed Brown         switch (user->testNum) {
708c4762a1bSJed Brown         case 0: {
709c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 2};
710c4762a1bSJed Brown           PetscInt tetPoints_p2[3] = {0, 1, 2};
711c4762a1bSJed Brown 
712c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 3, &points);CHKERRQ(ierr);
713c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  tetSizes_p2, 2);CHKERRQ(ierr);
714c4762a1bSJed Brown           ierr = PetscArraycpy(points, tetPoints_p2, 3);CHKERRQ(ierr);break;}
715c4762a1bSJed Brown         case 1: {
716c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {3, 4};
717c4762a1bSJed Brown           PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6};
718c4762a1bSJed Brown 
719c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 7, &points);CHKERRQ(ierr);
720c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  tetSizes_p2, 2);CHKERRQ(ierr);
721c4762a1bSJed Brown           ierr = PetscArraycpy(points, tetPoints_p2, 7);CHKERRQ(ierr);break;}
722c4762a1bSJed Brown         default:
723c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for tetrahedral hybrid mesh on 2 procs", user->testNum);
724c4762a1bSJed Brown         }
725c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) {
726c4762a1bSJed Brown         switch (user->testNum) {
727c4762a1bSJed Brown         case 0: {
728c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 1};
729c4762a1bSJed Brown           PetscInt hexPoints_p2[2] = {0, 1};
730c4762a1bSJed Brown 
731c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
732c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  hexSizes_p2, 2);CHKERRQ(ierr);
733c4762a1bSJed Brown           ierr = PetscArraycpy(points, hexPoints_p2, 2);CHKERRQ(ierr);break;}
734c4762a1bSJed Brown         default:
735c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral mesh on 2 procs", user->testNum);
736c4762a1bSJed Brown         }
737c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) {
738c4762a1bSJed Brown         switch (user->testNum) {
739c4762a1bSJed Brown         case 0: {
740c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 1};
741c4762a1bSJed Brown           PetscInt hexPoints_p2[2] = {0, 1};
742c4762a1bSJed Brown 
743c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 2, &points);CHKERRQ(ierr);
744c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  hexSizes_p2, 2);CHKERRQ(ierr);
745c4762a1bSJed Brown           ierr = PetscArraycpy(points, hexPoints_p2, 2);CHKERRQ(ierr);break;}
746c4762a1bSJed Brown         case 1: {
747c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {5, 4};
748c4762a1bSJed Brown           PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6};
749c4762a1bSJed Brown 
750c4762a1bSJed Brown           ierr = PetscMalloc2(2, &sizes, 9, &points);CHKERRQ(ierr);
751c4762a1bSJed Brown           ierr = PetscArraycpy(sizes,  hexSizes_p2, 2);CHKERRQ(ierr);
752c4762a1bSJed Brown           ierr = PetscArraycpy(points, hexPoints_p2, 9);CHKERRQ(ierr);break;}
753c4762a1bSJed Brown         default:
754c4762a1bSJed Brown           SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %d for hexahedral hybrid mesh on 2 procs", user->testNum);
755c4762a1bSJed Brown         }
756c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
757c4762a1bSJed Brown     }
758c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm, &part);CHKERRQ(ierr);
759c4762a1bSJed Brown     ierr = PetscPartitionerSetType(part, PETSCPARTITIONERSHELL);CHKERRQ(ierr);
760c4762a1bSJed Brown     ierr = PetscPartitionerShellSetPartition(part, size, sizes, points);CHKERRQ(ierr);
761c4762a1bSJed Brown     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
762c4762a1bSJed Brown   } else {
763c4762a1bSJed Brown     PetscPartitioner part;
764c4762a1bSJed Brown 
765c4762a1bSJed Brown     ierr = DMPlexGetPartitioner(*dm,&part);CHKERRQ(ierr);
766c4762a1bSJed Brown     ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
767c4762a1bSJed Brown   }
768c4762a1bSJed Brown   {
769c4762a1bSJed Brown     DM pdm = NULL;
770c4762a1bSJed Brown 
771c4762a1bSJed Brown     ierr = DMPlexDistribute(*dm, 0, NULL, &pdm);CHKERRQ(ierr);
772c4762a1bSJed Brown     if (pdm) {
773c4762a1bSJed Brown       ierr = DMViewFromOptions(pdm, NULL, "-dm_view");CHKERRQ(ierr);
774c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
775c4762a1bSJed Brown       *dm  = pdm;
776c4762a1bSJed Brown     }
777c4762a1bSJed Brown   }
778c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
779c4762a1bSJed Brown   if (user->uninterpolate || user->reinterpolate) {
780c4762a1bSJed Brown     DM udm = NULL;
781c4762a1bSJed Brown 
782c4762a1bSJed Brown     ierr = DMPlexUninterpolate(*dm, &udm);CHKERRQ(ierr);
783c4762a1bSJed Brown     ierr = DMPlexCopyCoordinates(*dm, udm);CHKERRQ(ierr);
784c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
785c4762a1bSJed Brown     *dm  = udm;
786c4762a1bSJed Brown   }
787c4762a1bSJed Brown   if (user->reinterpolate) {
788c4762a1bSJed Brown     DM idm = NULL;
789c4762a1bSJed Brown 
790c4762a1bSJed Brown     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
791c4762a1bSJed Brown     ierr = DMPlexCopyCoordinates(*dm, idm);CHKERRQ(ierr);
792c4762a1bSJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
793c4762a1bSJed Brown     *dm  = idm;
794c4762a1bSJed Brown   }
795c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Hybrid Mesh");CHKERRQ(ierr);
796c4762a1bSJed Brown   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
797c4762a1bSJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject) *dm, "hyb_");CHKERRQ(ierr);
798c4762a1bSJed Brown   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
799c4762a1bSJed Brown   PetscFunctionReturn(0);
800c4762a1bSJed Brown }
801c4762a1bSJed Brown 
802c4762a1bSJed Brown int main(int argc, char **argv)
803c4762a1bSJed Brown {
804c4762a1bSJed Brown   DM             dm;
805c4762a1bSJed Brown   AppCtx         user;                 /* user-defined work context */
806c4762a1bSJed Brown   PetscErrorCode ierr;
807c4762a1bSJed Brown 
808c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
809c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
810c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
811c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
812c4762a1bSJed Brown   ierr = PetscFinalize();
813c4762a1bSJed Brown   return ierr;
814c4762a1bSJed Brown }
815c4762a1bSJed Brown 
816c4762a1bSJed Brown /*TEST
817c4762a1bSJed Brown 
818c4762a1bSJed Brown   # 1D Simplex 29-31
819c4762a1bSJed Brown   testset:
820*54fcfd0cSMatthew G. Knepley     args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
821c4762a1bSJed Brown     test:
822c4762a1bSJed Brown       suffix: 29
823c4762a1bSJed Brown     test:
824c4762a1bSJed Brown       suffix: 30
825c4762a1bSJed Brown       args: -dm_refine 1
826c4762a1bSJed Brown     test:
827c4762a1bSJed Brown       suffix: 31
828c4762a1bSJed Brown       args: -dm_refine 5
829c4762a1bSJed Brown 
830c4762a1bSJed Brown   # 2D Simplex 0-3
831c4762a1bSJed Brown   testset:
832*54fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
833c4762a1bSJed Brown     test:
834c4762a1bSJed Brown       suffix: 0
835c4762a1bSJed Brown       args: -hyb_dm_plex_check_faces
836c4762a1bSJed Brown     test:
837c4762a1bSJed Brown       suffix: 1
838c4762a1bSJed Brown       args: -dm_refine 1 -hyb_dm_plex_check_faces
839c4762a1bSJed Brown     test:
840c4762a1bSJed Brown       suffix: 2
841c4762a1bSJed Brown       nsize: 2
842c4762a1bSJed Brown       args: -hyb_dm_plex_check_faces
843c4762a1bSJed Brown     test:
844c4762a1bSJed Brown       suffix: 3
845c4762a1bSJed Brown       nsize: 2
846c4762a1bSJed Brown       args: -dm_refine 1 -hyb_dm_plex_check_faces
847c4762a1bSJed Brown     test:
848c4762a1bSJed Brown       suffix: 32
849c4762a1bSJed Brown       args: -dm_refine 1 -uninterpolate
850c4762a1bSJed Brown     test:
851c4762a1bSJed Brown       suffix: 33
852c4762a1bSJed Brown       nsize: 2
853c4762a1bSJed Brown       args: -dm_refine 1 -uninterpolate
854c4762a1bSJed Brown     test:
855c4762a1bSJed Brown       suffix: 34
856c4762a1bSJed Brown       nsize: 2
857c4762a1bSJed Brown       args: -dm_refine 3 -uninterpolate
858c4762a1bSJed Brown 
859c4762a1bSJed Brown   # 2D Hybrid Simplex 4-7
860c4762a1bSJed Brown   testset:
861*54fcfd0cSMatthew G. Knepley     args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
862c4762a1bSJed Brown     test:
863c4762a1bSJed Brown       suffix: 4
864c4762a1bSJed Brown     test:
865c4762a1bSJed Brown       suffix: 5
866c4762a1bSJed Brown       args: -dm_refine 1
867c4762a1bSJed Brown     test:
868c4762a1bSJed Brown       suffix: 6
869c4762a1bSJed Brown       nsize: 2
870c4762a1bSJed Brown     test:
871c4762a1bSJed Brown       suffix: 7
872c4762a1bSJed Brown       nsize: 2
873c4762a1bSJed Brown       args: -dm_refine 1
874c4762a1bSJed Brown     test:
875c4762a1bSJed Brown       suffix: 24
876c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
877c4762a1bSJed Brown 
878c4762a1bSJed Brown   # 2D Quad 12-13
879c4762a1bSJed Brown   testset:
880*54fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
881c4762a1bSJed Brown     test:
882c4762a1bSJed Brown       suffix: 12
883c4762a1bSJed Brown       args: -dm_refine 1
884c4762a1bSJed Brown     test:
885c4762a1bSJed Brown       suffix: 13
886c4762a1bSJed Brown       nsize: 2
887c4762a1bSJed Brown       args: -dm_refine 1
888c4762a1bSJed Brown 
889c4762a1bSJed Brown   # 2D Hybrid Quad 27-28
890c4762a1bSJed Brown   testset:
891*54fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
892c4762a1bSJed Brown     test:
893c4762a1bSJed Brown       suffix: 27
894c4762a1bSJed Brown       args: -dm_refine 1
895c4762a1bSJed Brown     test:
896c4762a1bSJed Brown       suffix: 28
897c4762a1bSJed Brown       nsize: 2
898c4762a1bSJed Brown       args: -dm_refine 1
899c4762a1bSJed Brown 
900c4762a1bSJed Brown   # 3D Simplex 8-11
901c4762a1bSJed Brown   testset:
902*54fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
903c4762a1bSJed Brown     test:
904c4762a1bSJed Brown       suffix: 8
905c4762a1bSJed Brown       args: -dm_refine 1
906c4762a1bSJed Brown     test:
907c4762a1bSJed Brown       suffix: 9
908c4762a1bSJed Brown       nsize: 2
909c4762a1bSJed Brown       args: -dm_refine 1
910c4762a1bSJed Brown     test:
911c4762a1bSJed Brown       suffix: 10
912c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
913c4762a1bSJed Brown     test:
914c4762a1bSJed Brown       suffix: 11
915c4762a1bSJed Brown       nsize: 2
916c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
917c4762a1bSJed Brown     test:
918c4762a1bSJed Brown       suffix: 25
919c4762a1bSJed Brown       args: -test_num 2 -dm_refine 2
920c4762a1bSJed Brown 
921c4762a1bSJed Brown   # 3D Hybrid Simplex 16-19
922c4762a1bSJed Brown   testset:
923*54fcfd0cSMatthew G. Knepley     args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
924c4762a1bSJed Brown     test:
925c4762a1bSJed Brown       suffix: 16
926c4762a1bSJed Brown       args: -dm_refine 1
927c4762a1bSJed Brown     test:
928c4762a1bSJed Brown       suffix: 17
929c4762a1bSJed Brown       nsize: 2
930c4762a1bSJed Brown       args: -dm_refine 1
931c4762a1bSJed Brown     test:
932c4762a1bSJed Brown       suffix: 18
933c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
934c4762a1bSJed Brown     test:
935c4762a1bSJed Brown       suffix: 19
936c4762a1bSJed Brown       nsize: 2
937c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
938c4762a1bSJed Brown 
939c4762a1bSJed Brown   # 3D Hex 14-15
940c4762a1bSJed Brown   testset:
941*54fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
942c4762a1bSJed Brown     test:
943c4762a1bSJed Brown       suffix: 14
944c4762a1bSJed Brown       args: -dm_refine 1
945c4762a1bSJed Brown     test:
946c4762a1bSJed Brown       suffix: 15
947c4762a1bSJed Brown       nsize: 2
948c4762a1bSJed Brown      args: -dm_refine 1
949c4762a1bSJed Brown     test:
950c4762a1bSJed Brown       suffix: 26
951c4762a1bSJed Brown       args: -test_num 1 -dm_refine 2
952c4762a1bSJed Brown 
953c4762a1bSJed Brown   # 3D Hybrid Hex 20-23
954c4762a1bSJed Brown   testset:
955*54fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
956c4762a1bSJed Brown     test:
957c4762a1bSJed Brown       suffix: 20
958c4762a1bSJed Brown       args: -dm_refine 1
959c4762a1bSJed Brown     test:
960c4762a1bSJed Brown       suffix: 21
961c4762a1bSJed Brown       nsize: 2
962c4762a1bSJed Brown       args: -dm_refine 1
963c4762a1bSJed Brown     test:
964c4762a1bSJed Brown       suffix: 22
965c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
966c4762a1bSJed Brown     test:
967c4762a1bSJed Brown       suffix: 23
968c4762a1bSJed Brown       nsize: 2
969c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
970c4762a1bSJed Brown 
971c4762a1bSJed Brown   # Hybrid interpolation
972*54fcfd0cSMatthew G. Knepley   #   TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells
973c4762a1bSJed Brown   testset:
974c4762a1bSJed Brown     nsize: 2
975*54fcfd0cSMatthew G. Knepley     args: -test_partition 0 -petscpartitioner_type simple -dm_view -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
976c4762a1bSJed Brown     test:
977c4762a1bSJed Brown       suffix: hybint_2d_0
978c4762a1bSJed Brown       args: -dim 2 -dm_refine 2
979c4762a1bSJed Brown     test:
980c4762a1bSJed Brown       suffix: hybint_2d_1
981c4762a1bSJed Brown       args: -dim 2 -dm_refine 2 -test_num 1
982c4762a1bSJed Brown     test:
983c4762a1bSJed Brown       suffix: hybint_3d_0
984c4762a1bSJed Brown       args: -dim 3 -dm_refine 1
985c4762a1bSJed Brown     test:
986c4762a1bSJed Brown       suffix: hybint_3d_1
987c4762a1bSJed Brown       args: -dim 3 -dm_refine 1 -test_num 1
988c4762a1bSJed Brown 
989c4762a1bSJed Brown TEST*/
990