xref: /petsc/src/dm/impls/plex/tests/ex4.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
154fcfd0cSMatthew 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 
ProcessOptions(MPI_Comm comm,AppCtx * options)16d71ae5a4SJacob Faibussowitsch PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
17d71ae5a4SJacob Faibussowitsch {
18c4762a1bSJed Brown   PetscFunctionBegin;
19c4762a1bSJed Brown   options->debug         = 0;
20c4762a1bSJed Brown   options->dim           = 2;
21c4762a1bSJed Brown   options->cellHybrid    = PETSC_TRUE;
22c4762a1bSJed Brown   options->cellSimplex   = PETSC_TRUE;
23c4762a1bSJed Brown   options->testPartition = PETSC_TRUE;
24c4762a1bSJed Brown   options->testNum       = 0;
25c4762a1bSJed Brown   options->uninterpolate = PETSC_FALSE;
26c4762a1bSJed Brown   options->reinterpolate = PETSC_FALSE;
27c4762a1bSJed Brown 
28d0609cedSBarry Smith   PetscOptionsBegin(comm, "", "Meshing Problem Options", "DMPLEX");
299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-debug", "The debugging level", "ex4.c", options->debug, &options->debug, NULL, 0));
309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex4.c", options->dim, &options->dim, NULL, 1, 3));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_hybrid", "Use a hybrid mesh", "ex4.c", options->cellHybrid, &options->cellHybrid, NULL));
329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-cell_simplex", "Use simplices if true, otherwise hexes", "ex4.c", options->cellSimplex, &options->cellSimplex, NULL));
339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-test_partition", "Use a fixed partition for testing", "ex4.c", options->testPartition, &options->testPartition, NULL));
349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-test_num", "The particular mesh to test", "ex4.c", options->testNum, &options->testNum, NULL, 0));
359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-uninterpolate", "Uninterpolate the mesh at the end", "ex4.c", options->uninterpolate, &options->uninterpolate, NULL));
369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-reinterpolate", "Reinterpolate the mesh at the end", "ex4.c", options->reinterpolate, &options->reinterpolate, NULL));
37d0609cedSBarry Smith   PetscOptionsEnd();
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown /* Two segments
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   2-------0-------3-------1-------4
44c4762a1bSJed Brown 
45c4762a1bSJed Brown become
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   4---0---7---1---5---2---8---3---6
48c4762a1bSJed Brown 
49c4762a1bSJed Brown */
CreateSimplex_1D(MPI_Comm comm,DM * dm)50d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplex_1D(MPI_Comm comm, DM *dm)
51d71ae5a4SJacob Faibussowitsch {
52c4762a1bSJed Brown   PetscInt    depth = 1;
53c4762a1bSJed Brown   PetscMPIInt rank;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
57dd400576SPatrick Sanan   if (rank == 0) {
58c4762a1bSJed Brown     PetscInt    numPoints[2]         = {3, 2};
59c4762a1bSJed Brown     PetscInt    coneSize[5]          = {2, 2, 0, 0, 0};
60c4762a1bSJed Brown     PetscInt    cones[4]             = {2, 3, 3, 4};
61c4762a1bSJed Brown     PetscInt    coneOrientations[16] = {0, 0, 0, 0};
62c4762a1bSJed Brown     PetscScalar vertexCoords[3]      = {-1.0, 0.0, 1.0};
63c4762a1bSJed Brown 
649566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
65c4762a1bSJed Brown   } else {
66c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
67c4762a1bSJed Brown 
689566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
69c4762a1bSJed Brown   }
703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71c4762a1bSJed Brown }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown /* Two triangles
74c4762a1bSJed Brown         4
75c4762a1bSJed Brown       / | \
76c4762a1bSJed Brown      8  |  10
77c4762a1bSJed Brown     /   |   \
78c4762a1bSJed Brown    2  0 7  1 5
79c4762a1bSJed Brown     \   |   /
80c4762a1bSJed Brown      6  |  9
81c4762a1bSJed Brown       \ | /
82c4762a1bSJed Brown         3
83c4762a1bSJed Brown 
84c4762a1bSJed Brown Becomes
85c4762a1bSJed Brown            10
86c4762a1bSJed Brown           / | \
87c4762a1bSJed Brown         21  |  26
88c4762a1bSJed Brown         /   |   \
89c4762a1bSJed Brown       14 2 20 4  16
90c4762a1bSJed Brown       /|\   |   /|\
91c4762a1bSJed Brown     22 | 28 | 32 | 25
92c4762a1bSJed Brown     /  |  \ | /  | 6\
93c4762a1bSJed Brown    8  29 3 13  7 31  11
94c4762a1bSJed Brown     \0 |  / | \  |  /
95c4762a1bSJed Brown     17 | 27 | 30 | 24
96c4762a1bSJed Brown       \|/   |   \|/
97c4762a1bSJed Brown       12 1 19 5  15
98c4762a1bSJed Brown         \   |   /
99c4762a1bSJed Brown         18  |  23
100c4762a1bSJed Brown           \ | /
101c4762a1bSJed Brown             9
102c4762a1bSJed Brown */
CreateSimplex_2D(MPI_Comm comm,DM * dm)103d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplex_2D(MPI_Comm comm, DM *dm)
104d71ae5a4SJacob Faibussowitsch {
105c4762a1bSJed Brown   PetscInt    depth = 2;
106c4762a1bSJed Brown   PetscMPIInt rank;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   PetscFunctionBegin;
1099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
110dd400576SPatrick Sanan   if (rank == 0) {
111c4762a1bSJed Brown     PetscInt    numPoints[3]         = {4, 5, 2};
112c4762a1bSJed Brown     PetscInt    coneSize[11]         = {3, 3, 0, 0, 0, 0, 2, 2, 2, 2, 2};
113c4762a1bSJed Brown     PetscInt    cones[16]            = {6, 7, 8, 7, 9, 10, 2, 3, 3, 4, 4, 2, 3, 5, 5, 4};
114b5a892a1SMatthew G. Knepley     PetscInt    coneOrientations[16] = {0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
115c4762a1bSJed Brown     PetscScalar vertexCoords[8]      = {-0.5, 0.0, 0.0, -0.5, 0.0, 0.5, 0.5, 0.0};
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
118c4762a1bSJed Brown   } else {
119c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
122c4762a1bSJed Brown   }
1233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /* Two triangles separated by a zero-volume cell with 4 vertices/2 edges
127c4762a1bSJed Brown         5--16--8
128c4762a1bSJed Brown       / |      | \
129c4762a1bSJed Brown     11  |      |  12
130c4762a1bSJed Brown     /   |      |   \
131c4762a1bSJed Brown    3  0 10  2 14 1  6
132c4762a1bSJed Brown     \   |      |   /
133c4762a1bSJed Brown      9  |      |  13
134c4762a1bSJed Brown       \ |      | /
135c4762a1bSJed Brown         4--15--7
136c4762a1bSJed Brown */
CreateSimplexHybrid_2D(MPI_Comm comm,PetscInt testNum,DM * dm)137d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplexHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
138d71ae5a4SJacob Faibussowitsch {
139c4762a1bSJed Brown   DM          idm, hdm = NULL;
140c4762a1bSJed Brown   DMLabel     faultLabel, hybridLabel;
141c4762a1bSJed Brown   PetscInt    p;
142c4762a1bSJed Brown   PetscMPIInt rank;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
146dd400576SPatrick Sanan   if (rank == 0) {
147c4762a1bSJed Brown     switch (testNum) {
1489371c9d4SSatish Balay     case 0: {
149c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 2};
150c4762a1bSJed Brown       PetscInt    coneSize[6]         = {3, 3, 0, 0, 0, 0};
151c4762a1bSJed Brown       PetscInt    cones[6]            = {2, 3, 4, 5, 4, 3};
152c4762a1bSJed Brown       PetscInt    coneOrientations[6] = {0, 0, 0, 0, 0, 0};
153c4762a1bSJed Brown       PetscScalar vertexCoords[8]     = {-1.0, -0.5, 0.0, -0.5, 0.0, 0.5, 1.0, 0.5};
154c4762a1bSJed Brown       PetscInt    faultPoints[2]      = {3, 4};
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
1579566063dSJacob Faibussowitsch       for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
1589371c9d4SSatish Balay     } break;
1599371c9d4SSatish Balay     case 1: {
160c4762a1bSJed Brown       PetscInt    numPoints[2]         = {5, 4};
161c4762a1bSJed Brown       PetscInt    coneSize[9]          = {3, 3, 3, 3, 0, 0, 0, 0, 0};
162c4762a1bSJed Brown       PetscInt    cones[12]            = {4, 5, 6, 6, 7, 4, 6, 5, 8, 6, 8, 7};
163c4762a1bSJed Brown       PetscInt    coneOrientations[12] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
164c4762a1bSJed Brown       PetscScalar vertexCoords[10]     = {-1.0, 0.0, 0.0, -1.0, 0.0, 0.0, 0.0, 1.0, 1.0, 0.0};
165c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {5, 6, 7};
166c4762a1bSJed Brown 
1679566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
1689566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
1699371c9d4SSatish Balay     } break;
170d71ae5a4SJacob Faibussowitsch     default:
171d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
172c4762a1bSJed Brown     }
1739566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
1749566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
1759566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
1769566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
1779566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
1789566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
179caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
1809566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
1819566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
1829566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
183c4762a1bSJed Brown     *dm = hdm;
184c4762a1bSJed Brown   } else {
185c4762a1bSJed Brown     PetscInt numPoints[2] = {0, 0};
186c4762a1bSJed Brown 
1879566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
1889566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
1899566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
1909566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
1919566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
1929566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
193caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
1949566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
1959566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
196c4762a1bSJed Brown     *dm = hdm;
197c4762a1bSJed Brown   }
1983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
199c4762a1bSJed Brown }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown /* Two quadrilaterals
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   5----10-----4----14-----7
204c4762a1bSJed Brown   |           |           |
205c4762a1bSJed Brown   |           |           |
206c4762a1bSJed Brown   |           |           |
207c4762a1bSJed Brown  11     0     9     1     13
208c4762a1bSJed Brown   |           |           |
209c4762a1bSJed Brown   |           |           |
210c4762a1bSJed Brown   |           |           |
211c4762a1bSJed Brown   2-----8-----3----12-----6
212c4762a1bSJed Brown */
CreateTensorProduct_2D(MPI_Comm comm,PetscInt testNum,DM * dm)213d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateTensorProduct_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
214d71ae5a4SJacob Faibussowitsch {
215c4762a1bSJed Brown   PetscInt    depth = 2;
216c4762a1bSJed Brown   PetscMPIInt rank;
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   PetscFunctionBegin;
2199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
220dd400576SPatrick Sanan   if (rank == 0) {
221c4762a1bSJed Brown     PetscInt    numPoints[3]         = {6, 7, 2};
222c4762a1bSJed Brown     PetscInt    coneSize[15]         = {4, 4, 0, 0, 0, 0, 0, 0, 2, 2, 2, 2, 2, 2, 2};
223c4762a1bSJed 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};
224b5a892a1SMatthew G. Knepley     PetscInt    coneOrientations[22] = {0, 0, 0, 0, 0, 0, 0, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
225c4762a1bSJed 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};
226c4762a1bSJed Brown 
2279566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
228c4762a1bSJed Brown   } else {
229c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
232c4762a1bSJed Brown   }
2333ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
234c4762a1bSJed Brown }
235c4762a1bSJed Brown 
CreateTensorProductHybrid_2D(MPI_Comm comm,PetscInt testNum,DM * dm)236d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateTensorProductHybrid_2D(MPI_Comm comm, PetscInt testNum, DM *dm)
237d71ae5a4SJacob Faibussowitsch {
238c4762a1bSJed Brown   DM          idm, hdm = NULL;
239c4762a1bSJed Brown   DMLabel     faultLabel, hybridLabel;
240c4762a1bSJed Brown   PetscInt    p;
241c4762a1bSJed Brown   PetscMPIInt rank;
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   PetscFunctionBegin;
2449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
245dd400576SPatrick Sanan   if (rank == 0) {
246c4762a1bSJed Brown     PetscInt numPoints[2] = {6, 2};
247c4762a1bSJed Brown     PetscInt coneSize[8]  = {4, 4, 0, 0, 0, 0, 0, 0};
2489371c9d4SSatish Balay     PetscInt cones[8]     = {
2499371c9d4SSatish Balay       2, 3, 4, 5, 3, 6, 7, 4,
2509371c9d4SSatish Balay     };
251c4762a1bSJed Brown     PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
252c4762a1bSJed 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};
253c4762a1bSJed Brown     PetscInt    faultPoints[2]      = {3, 4};
254c4762a1bSJed Brown 
2559566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
2569566063dSJacob Faibussowitsch     for (p = 0; p < 2; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
2579566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
2589566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
2599566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
2609566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
2619566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
2629566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
263caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
2649566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
265c4762a1bSJed Brown   } else {
266c4762a1bSJed Brown     PetscInt numPoints[3] = {0, 0, 0};
267c4762a1bSJed Brown 
2689566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
2699566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
2709566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
2719566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
2729566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
2739566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
274caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
275c4762a1bSJed Brown   }
2769566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&idm));
2779566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
278c4762a1bSJed Brown   *dm = hdm;
2793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
280c4762a1bSJed Brown }
281c4762a1bSJed Brown 
282c4762a1bSJed Brown /* Two tetrahedrons
283c4762a1bSJed Brown 
284c4762a1bSJed Brown  cell   5          5______    cell
285c4762a1bSJed Brown  0    / | \        |\      \     1
286c4762a1bSJed Brown     17  |  18      | 18 13  21
287c4762a1bSJed Brown     /8 19 10\     19  \      \
288c4762a1bSJed Brown    2-14-|----4     |   4--22--6
289c4762a1bSJed Brown     \ 9 | 7 /      |10 /      /
290c4762a1bSJed Brown     16  |  15      | 15  12 20
291c4762a1bSJed Brown       \ | /        |/      /
292c4762a1bSJed Brown         3          3------
293c4762a1bSJed Brown */
CreateSimplex_3D(MPI_Comm comm,PetscInt testNum,DM * dm)294d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplex_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
295d71ae5a4SJacob Faibussowitsch {
296c4762a1bSJed Brown   DM          idm;
297c4762a1bSJed Brown   PetscInt    depth = 3;
298c4762a1bSJed Brown   PetscMPIInt rank;
299c4762a1bSJed Brown 
300c4762a1bSJed Brown   PetscFunctionBegin;
3019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
302dd400576SPatrick Sanan   if (rank == 0) {
303c4762a1bSJed Brown     switch (testNum) {
3049371c9d4SSatish Balay     case 0: {
305c4762a1bSJed Brown       PetscInt    numPoints[4]         = {5, 9, 7, 2};
306c4762a1bSJed 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};
307c4762a1bSJed 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};
308b5a892a1SMatthew G. Knepley       PetscInt    coneOrientations[47] = {0, 0, 0, 0, -2, 0, 0, 0, 0, 0, 0, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, -1, -1, 0, -1, -1, -1, -1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
309c4762a1bSJed 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};
310c4762a1bSJed Brown 
3119566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3129371c9d4SSatish Balay     } break;
3139371c9d4SSatish Balay     case 1: {
314c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
315c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
316c4762a1bSJed Brown       PetscInt    cones[8]            = {4, 3, 5, 2, 5, 3, 4, 6};
317c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
318c4762a1bSJed 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};
319c4762a1bSJed Brown 
320c4762a1bSJed Brown       depth = 1;
3219566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3229566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
3239566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
3249566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
325c4762a1bSJed Brown       *dm = idm;
3269371c9d4SSatish Balay     } break;
3279371c9d4SSatish Balay     case 2: {
328c4762a1bSJed Brown       PetscInt    numPoints[2]        = {4, 1};
329c4762a1bSJed Brown       PetscInt    coneSize[5]         = {4, 0, 0, 0, 0};
330c4762a1bSJed Brown       PetscInt    cones[4]            = {2, 3, 4, 1};
331c4762a1bSJed Brown       PetscInt    coneOrientations[4] = {0, 0, 0, 0};
332c4762a1bSJed 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};
333c4762a1bSJed Brown 
334c4762a1bSJed Brown       depth = 1;
3359566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3369566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
3379566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
3389566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
339c4762a1bSJed Brown       *dm = idm;
3409371c9d4SSatish Balay     } break;
341d71ae5a4SJacob Faibussowitsch     default:
342d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
343c4762a1bSJed Brown     }
344c4762a1bSJed Brown   } else {
345c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
346c4762a1bSJed Brown 
3479566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, depth, numPoints, NULL, NULL, NULL, NULL));
348c4762a1bSJed Brown     switch (testNum) {
349c4762a1bSJed Brown     case 1:
3509566063dSJacob Faibussowitsch       PetscCall(DMPlexInterpolate(*dm, &idm));
3519566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
3529566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
353c4762a1bSJed Brown       *dm = idm;
354c4762a1bSJed Brown       break;
355c4762a1bSJed Brown     }
356c4762a1bSJed Brown   }
3573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
358c4762a1bSJed Brown }
359c4762a1bSJed Brown 
360c4762a1bSJed Brown /* Two tetrahedrons separated by a zero-volume cell with 6 vertices
361c4762a1bSJed Brown 
362c4762a1bSJed Brown  cell   6 ___33___10______    cell
363c4762a1bSJed Brown  0    / | \        |\      \     1
364c4762a1bSJed Brown     21  |  23      | 29     27
365c4762a1bSJed Brown     /12 24 14\    30  \      \
366c4762a1bSJed Brown    3-20-|----5--32-|---9--26--7
367c4762a1bSJed Brown     \ 13| 11/      |18 /      /
368c4762a1bSJed Brown     19  |  22      | 28     25
369c4762a1bSJed Brown       \ | /        |/      /
370c4762a1bSJed Brown         4----31----8------
371c4762a1bSJed Brown          cell 2
372c4762a1bSJed Brown */
CreateSimplexHybrid_3D(MPI_Comm comm,PetscInt testNum,DM * dm)373d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateSimplexHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
374d71ae5a4SJacob Faibussowitsch {
375c4762a1bSJed Brown   DM          idm, hdm = NULL;
376c4762a1bSJed Brown   DMLabel     faultLabel, hybridLabel;
377c4762a1bSJed Brown   PetscInt    p;
378c4762a1bSJed Brown   PetscMPIInt rank;
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
382dd400576SPatrick Sanan   if (rank == 0) {
383c4762a1bSJed Brown     switch (testNum) {
3849371c9d4SSatish Balay     case 0: {
385c4762a1bSJed Brown       PetscInt    numPoints[2]        = {5, 2};
386c4762a1bSJed Brown       PetscInt    coneSize[7]         = {4, 4, 0, 0, 0, 0, 0};
387c4762a1bSJed Brown       PetscInt    cones[8]            = {4, 3, 5, 2, 5, 3, 4, 6};
388c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
389c4762a1bSJed 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};
390c4762a1bSJed Brown       PetscInt    faultPoints[3]      = {3, 4, 5};
391c4762a1bSJed Brown 
3929566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
3939566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
3949371c9d4SSatish Balay     } break;
3959371c9d4SSatish Balay     case 1: {
396c4762a1bSJed Brown       /* Tets 0,3,5 and 1,2,4 */
397c4762a1bSJed Brown       PetscInt    numPoints[2]         = {9, 6};
398c4762a1bSJed Brown       PetscInt    coneSize[15]         = {4, 4, 4, 4, 4, 4, 0, 0, 0, 0, 0, 0, 0, 0, 0};
3999371c9d4SSatish Balay       PetscInt    cones[24]            = {7, 8, 9, 6, 11, 13, 9, 14, 10, 11, 13, 9, 9, 10, 11, 7, 9, 14, 13, 12, 7, 8, 11, 9};
4009371c9d4SSatish Balay       PetscInt    coneOrientations[24] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4019371c9d4SSatish Balay       PetscScalar vertexCoords[27]     = {-2.0, -1.0, 0.0, -2.0, 0.0, 0.0, -2.0, 0.0, 1.0, 0.0, -1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 2.0, -1.0, 0.0, 2.0, 0.0, 0.0, 2.0, 0.0, 1.0};
402c4762a1bSJed Brown       PetscInt    faultPoints[3]       = {9, 10, 11};
403c4762a1bSJed Brown 
4049566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4059566063dSJacob Faibussowitsch       for (p = 0; p < 3; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
4069371c9d4SSatish Balay     } break;
407d71ae5a4SJacob Faibussowitsch     default:
408d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
409c4762a1bSJed Brown     }
4109566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
4119566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
4129566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
4139566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
4149566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
4159566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
416caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, &hybridLabel, NULL, NULL, &hdm));
4179566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&hybridLabel));
4189566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
4199566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
420c4762a1bSJed Brown     *dm = hdm;
421c4762a1bSJed Brown   } else {
422c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
423c4762a1bSJed Brown 
4249566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
4259566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
4269566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
4279566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
4289566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
4299566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
430caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
4319566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
4329566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
433c4762a1bSJed Brown     *dm = hdm;
434c4762a1bSJed Brown   }
4353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
436c4762a1bSJed Brown }
437c4762a1bSJed Brown 
CreateTensorProduct_3D(MPI_Comm comm,PetscInt testNum,DM * dm)438d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateTensorProduct_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
439d71ae5a4SJacob Faibussowitsch {
440c4762a1bSJed Brown   DM          idm;
441c4762a1bSJed Brown   PetscMPIInt rank;
442c4762a1bSJed Brown 
443c4762a1bSJed Brown   PetscFunctionBegin;
4449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
445dd400576SPatrick Sanan   if (rank == 0) {
446c4762a1bSJed Brown     switch (testNum) {
4479371c9d4SSatish Balay     case 0: {
448c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
449c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
450c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8};
451c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4529371c9d4SSatish Balay       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, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5};
453c4762a1bSJed Brown 
4549566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4559371c9d4SSatish Balay     } break;
4569371c9d4SSatish Balay     case 1: {
457c4762a1bSJed Brown       PetscInt    numPoints[2]        = {8, 1};
458c4762a1bSJed Brown       PetscInt    coneSize[9]         = {8, 0, 0, 0, 0, 0, 0, 0, 0};
459c4762a1bSJed Brown       PetscInt    cones[8]            = {1, 2, 3, 4, 5, 6, 7, 8};
460c4762a1bSJed Brown       PetscInt    coneOrientations[8] = {0, 0, 0, 0, 0, 0, 0, 0};
4619371c9d4SSatish Balay       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, -1.0, -1.0, 1.0, 1.0, -1.0, 1.0, 1.0, 1.0, 1.0, -1.0, 1.0, 1.0};
462c4762a1bSJed Brown 
4639566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
4649371c9d4SSatish Balay     } break;
465d71ae5a4SJacob Faibussowitsch     default:
466d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
467c4762a1bSJed Brown     }
468c4762a1bSJed Brown   } else {
469c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
470c4762a1bSJed Brown 
4719566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
472c4762a1bSJed Brown   }
4739566063dSJacob Faibussowitsch   PetscCall(DMPlexInterpolate(*dm, &idm));
4749566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(idm, NULL, "-in_dm_view"));
4759566063dSJacob Faibussowitsch   PetscCall(DMDestroy(dm));
476c4762a1bSJed Brown   *dm = idm;
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
478c4762a1bSJed Brown }
479c4762a1bSJed Brown 
CreateTensorProductHybrid_3D(MPI_Comm comm,PetscInt testNum,DM * dm)480d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateTensorProductHybrid_3D(MPI_Comm comm, PetscInt testNum, DM *dm)
481d71ae5a4SJacob Faibussowitsch {
482c4762a1bSJed Brown   DM          idm, hdm = NULL;
483c4762a1bSJed Brown   DMLabel     faultLabel;
484c4762a1bSJed Brown   PetscInt    p;
485c4762a1bSJed Brown   PetscMPIInt rank;
486c4762a1bSJed Brown 
487c4762a1bSJed Brown   PetscFunctionBegin;
4889566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
489dd400576SPatrick Sanan   if (rank == 0) {
490c4762a1bSJed Brown     switch (testNum) {
4919371c9d4SSatish Balay     case 0: {
492c4762a1bSJed Brown       PetscInt    numPoints[2]         = {12, 2};
493c4762a1bSJed Brown       PetscInt    coneSize[14]         = {8, 8, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
494c4762a1bSJed Brown       PetscInt    cones[16]            = {2, 3, 4, 5, 6, 7, 8, 9, 5, 4, 10, 11, 7, 12, 13, 8};
495c4762a1bSJed Brown       PetscInt    coneOrientations[16] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
4969371c9d4SSatish Balay       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, -1.0, -0.5, 0.5, 0.0, -0.5, 0.5, 0.0, 0.5, 0.5, -1.0, 0.5, 0.5, 1.0, 0.5, -0.5, 1.0, -0.5, -0.5, 1.0, -0.5, 0.5, 1.0, 0.5, 0.5};
497c4762a1bSJed Brown       PetscInt    faultPoints[4]       = {2, 3, 5, 6};
498c4762a1bSJed Brown 
4999566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
5009566063dSJacob Faibussowitsch       for (p = 0; p < 4; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
5019371c9d4SSatish Balay     } break;
5029371c9d4SSatish Balay     case 1: {
503c4762a1bSJed Brown       PetscInt numPoints[2] = {30, 7};
504c4762a1bSJed 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};
5059371c9d4SSatish Balay       PetscInt cones[56] = {8, 21, 20, 7, 13, 12, 23, 24, 14, 15, 10, 9, 13, 8, 21, 24, 15, 16, 11, 10, 24, 21, 22, 25, 30, 29, 28, 21, 35, 24, 33, 34, 24, 21, 30, 35, 25, 36, 31, 22, 27, 20, 21, 28, 32, 33, 24, 23, 15, 24, 13, 14, 19, 18, 17, 26};
506c4762a1bSJed 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};
5079371c9d4SSatish Balay       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,  -2.0, -1.0, 0.0, -3.0, 0.0, 0.0, -2.0, 1.0, 0.0, -2.0, 2.0, 0.0,
5089371c9d4SSatish Balay                                        -2.0, -1.0, 2.0,  -3.0, 0.0,  2.0,  -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, 0.0,  2.0, 0.0, 0.0,  0.0, 2.0,
5099371c9d4SSatish Balay                                        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,  2.0,  -1.0, 0.0, 3.0,  0.0, 0.0, 2.0,  1.0, 0.0, 2.0,  2.0, 0.0};
510c4762a1bSJed Brown       PetscInt    faultPoints[6]    = {20, 21, 22, 23, 24, 25};
511c4762a1bSJed Brown 
5129566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords));
5139566063dSJacob Faibussowitsch       for (p = 0; p < 6; ++p) PetscCall(DMSetLabelValue(*dm, "fault", faultPoints[p], 1));
5149371c9d4SSatish Balay     } break;
515d71ae5a4SJacob Faibussowitsch     default:
516d71ae5a4SJacob Faibussowitsch       SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "No test mesh %" PetscInt_FMT, testNum);
517c4762a1bSJed Brown     }
5189566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
5199566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
5209566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
5219566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
5229566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
5239566063dSJacob Faibussowitsch     PetscCall(DMGetLabel(*dm, "fault", &faultLabel));
524caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, faultLabel, NULL, 0, NULL, NULL, NULL, &hdm));
5259566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
5269566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
527c4762a1bSJed Brown     *dm = hdm;
528c4762a1bSJed Brown   } else {
529c4762a1bSJed Brown     PetscInt numPoints[4] = {0, 0, 0, 0};
530c4762a1bSJed Brown 
5319566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFromDAG(*dm, 1, numPoints, NULL, NULL, NULL, NULL));
5329566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
5339566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)idm, "in_"));
5349566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeSetDefault(idm, PETSC_FALSE));
5359566063dSJacob Faibussowitsch     PetscCall(DMSetFromOptions(idm));
5369566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(idm, NULL, "-dm_view"));
537caf9e14dSMatthew G. Knepley     PetscCall(DMPlexCreateHybridMesh(idm, NULL, NULL, 0, NULL, NULL, NULL, &hdm));
5389566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&idm));
5399566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
540c4762a1bSJed Brown     *dm = hdm;
541c4762a1bSJed Brown   }
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
543c4762a1bSJed Brown }
544c4762a1bSJed Brown 
CreateMesh(MPI_Comm comm,AppCtx * user,DM * dm)545d71ae5a4SJacob Faibussowitsch PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
546d71ae5a4SJacob Faibussowitsch {
547c4762a1bSJed Brown   PetscInt    dim         = user->dim;
548c4762a1bSJed Brown   PetscBool   cellHybrid  = user->cellHybrid;
549c4762a1bSJed Brown   PetscBool   cellSimplex = user->cellSimplex;
550c4762a1bSJed Brown   PetscMPIInt rank, size;
551c4762a1bSJed Brown 
552c4762a1bSJed Brown   PetscFunctionBegin;
5539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5559566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
5569566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
5579566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
558c4762a1bSJed Brown   switch (dim) {
559c4762a1bSJed Brown   case 1:
56063a3b9bcSJacob Faibussowitsch     PetscCheck(!cellHybrid, comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make hybrid meshes for dimension %" PetscInt_FMT, dim);
5619566063dSJacob Faibussowitsch     PetscCall(CreateSimplex_1D(comm, dm));
562c4762a1bSJed Brown     break;
563c4762a1bSJed Brown   case 2:
564c4762a1bSJed Brown     if (cellSimplex) {
565c4762a1bSJed Brown       if (cellHybrid) {
5669566063dSJacob Faibussowitsch         PetscCall(CreateSimplexHybrid_2D(comm, user->testNum, dm));
567c4762a1bSJed Brown       } else {
5689566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_2D(comm, dm));
569c4762a1bSJed Brown       }
570c4762a1bSJed Brown     } else {
571c4762a1bSJed Brown       if (cellHybrid) {
5729566063dSJacob Faibussowitsch         PetscCall(CreateTensorProductHybrid_2D(comm, user->testNum, dm));
573c4762a1bSJed Brown       } else {
5749566063dSJacob Faibussowitsch         PetscCall(CreateTensorProduct_2D(comm, user->testNum, dm));
575c4762a1bSJed Brown       }
576c4762a1bSJed Brown     }
577c4762a1bSJed Brown     break;
578c4762a1bSJed Brown   case 3:
579c4762a1bSJed Brown     if (cellSimplex) {
580c4762a1bSJed Brown       if (cellHybrid) {
5819566063dSJacob Faibussowitsch         PetscCall(CreateSimplexHybrid_3D(comm, user->testNum, dm));
582c4762a1bSJed Brown       } else {
5839566063dSJacob Faibussowitsch         PetscCall(CreateSimplex_3D(comm, user->testNum, dm));
584c4762a1bSJed Brown       }
585c4762a1bSJed Brown     } else {
586c4762a1bSJed Brown       if (cellHybrid) {
5879566063dSJacob Faibussowitsch         PetscCall(CreateTensorProductHybrid_3D(comm, user->testNum, dm));
588c4762a1bSJed Brown       } else {
5899566063dSJacob Faibussowitsch         PetscCall(CreateTensorProduct_3D(comm, user->testNum, dm));
590c4762a1bSJed Brown       }
591c4762a1bSJed Brown     }
592c4762a1bSJed Brown     break;
593d71ae5a4SJacob Faibussowitsch   default:
594d71ae5a4SJacob Faibussowitsch     SETERRQ(comm, PETSC_ERR_ARG_OUTOFRANGE, "Cannot make meshes for dimension %" PetscInt_FMT, dim);
595c4762a1bSJed Brown   }
596c4762a1bSJed Brown   if (user->testPartition && size > 1) {
597c4762a1bSJed Brown     PetscPartitioner part;
598c4762a1bSJed Brown     PetscInt        *sizes  = NULL;
599c4762a1bSJed Brown     PetscInt        *points = NULL;
600c4762a1bSJed Brown 
601dd400576SPatrick Sanan     if (rank == 0) {
602c4762a1bSJed Brown       if (dim == 2 && cellSimplex && !cellHybrid && size == 2) {
603c4762a1bSJed Brown         switch (user->testNum) {
604c4762a1bSJed Brown         case 0: {
605c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 1};
606c4762a1bSJed Brown           PetscInt triPoints_p2[2] = {0, 1};
607c4762a1bSJed Brown 
6089566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
6099566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
6109371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, triPoints_p2, 2));
6119371c9d4SSatish Balay           break;
6129371c9d4SSatish Balay         }
613d71ae5a4SJacob Faibussowitsch         default:
614d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular mesh on 2 procs", user->testNum);
615c4762a1bSJed Brown         }
616c4762a1bSJed Brown       } else if (dim == 2 && cellSimplex && cellHybrid && size == 2) {
617c4762a1bSJed Brown         switch (user->testNum) {
618c4762a1bSJed Brown         case 0: {
619c4762a1bSJed Brown           PetscInt triSizes_p2[2]  = {1, 2};
620c4762a1bSJed Brown           PetscInt triPoints_p2[3] = {0, 1, 2};
621c4762a1bSJed Brown 
6229566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
6239566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, triSizes_p2, 2));
6249371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, triPoints_p2, 3));
6259371c9d4SSatish Balay           break;
6269371c9d4SSatish Balay         }
627d71ae5a4SJacob Faibussowitsch         default:
628d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for triangular hybrid mesh on 2 procs", user->testNum);
629c4762a1bSJed Brown         }
630c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && !cellHybrid && size == 2) {
631c4762a1bSJed Brown         switch (user->testNum) {
632c4762a1bSJed Brown         case 0: {
633c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 1};
634c4762a1bSJed Brown           PetscInt quadPoints_p2[2] = {0, 1};
635c4762a1bSJed Brown 
6369566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
6379566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
6389371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, quadPoints_p2, 2));
6399371c9d4SSatish Balay           break;
6409371c9d4SSatish Balay         }
641d71ae5a4SJacob Faibussowitsch         default:
642d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral mesh on 2 procs", user->testNum);
643c4762a1bSJed Brown         }
644c4762a1bSJed Brown       } else if (dim == 2 && !cellSimplex && cellHybrid && size == 2) {
645c4762a1bSJed Brown         switch (user->testNum) {
646c4762a1bSJed Brown         case 0: {
647c4762a1bSJed Brown           PetscInt quadSizes_p2[2]  = {1, 2};
648c4762a1bSJed Brown           PetscInt quadPoints_p2[3] = {0, 1, 2};
649c4762a1bSJed Brown 
6509566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
6519566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, quadSizes_p2, 2));
6529371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, quadPoints_p2, 3));
6539371c9d4SSatish Balay           break;
6549371c9d4SSatish Balay         }
655d71ae5a4SJacob Faibussowitsch         default:
656d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for quadrilateral hybrid mesh on 2 procs", user->testNum);
657c4762a1bSJed Brown         }
658c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && !cellHybrid && size == 2) {
659c4762a1bSJed Brown         switch (user->testNum) {
660c4762a1bSJed Brown         case 0: {
661c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 1};
662c4762a1bSJed Brown           PetscInt tetPoints_p2[2] = {0, 1};
663c4762a1bSJed Brown 
6649566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
6659566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
6669371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, tetPoints_p2, 2));
6679371c9d4SSatish Balay           break;
6689371c9d4SSatish Balay         }
669c4762a1bSJed Brown         case 1: {
670c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 1};
671c4762a1bSJed Brown           PetscInt tetPoints_p2[2] = {0, 1};
672c4762a1bSJed Brown 
6739566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
6749566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
6759371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, tetPoints_p2, 2));
6769371c9d4SSatish Balay           break;
6779371c9d4SSatish Balay         }
678d71ae5a4SJacob Faibussowitsch         default:
679d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral mesh on 2 procs", user->testNum);
680c4762a1bSJed Brown         }
681c4762a1bSJed Brown       } else if (dim == 3 && cellSimplex && cellHybrid && size == 2) {
682c4762a1bSJed Brown         switch (user->testNum) {
683c4762a1bSJed Brown         case 0: {
684c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {1, 2};
685c4762a1bSJed Brown           PetscInt tetPoints_p2[3] = {0, 1, 2};
686c4762a1bSJed Brown 
6879566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 3, &points));
6889566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
6899371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, tetPoints_p2, 3));
6909371c9d4SSatish Balay           break;
6919371c9d4SSatish Balay         }
692c4762a1bSJed Brown         case 1: {
693c4762a1bSJed Brown           PetscInt tetSizes_p2[2]  = {3, 4};
694c4762a1bSJed Brown           PetscInt tetPoints_p2[7] = {0, 3, 5, 1, 2, 4, 6};
695c4762a1bSJed Brown 
6969566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 7, &points));
6979566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, tetSizes_p2, 2));
6989371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, tetPoints_p2, 7));
6999371c9d4SSatish Balay           break;
7009371c9d4SSatish Balay         }
701d71ae5a4SJacob Faibussowitsch         default:
702d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for tetrahedral hybrid mesh on 2 procs", user->testNum);
703c4762a1bSJed Brown         }
704c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && !cellHybrid && size == 2) {
705c4762a1bSJed Brown         switch (user->testNum) {
706c4762a1bSJed Brown         case 0: {
707c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 1};
708c4762a1bSJed Brown           PetscInt hexPoints_p2[2] = {0, 1};
709c4762a1bSJed Brown 
7109566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
7119566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
7129371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, hexPoints_p2, 2));
7139371c9d4SSatish Balay           break;
7149371c9d4SSatish Balay         }
715d71ae5a4SJacob Faibussowitsch         default:
716d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral mesh on 2 procs", user->testNum);
717c4762a1bSJed Brown         }
718c4762a1bSJed Brown       } else if (dim == 3 && !cellSimplex && cellHybrid && size == 2) {
719c4762a1bSJed Brown         switch (user->testNum) {
720c4762a1bSJed Brown         case 0: {
721c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {1, 1};
722c4762a1bSJed Brown           PetscInt hexPoints_p2[2] = {0, 1};
723c4762a1bSJed Brown 
7249566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 2, &points));
7259566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
7269371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, hexPoints_p2, 2));
7279371c9d4SSatish Balay           break;
7289371c9d4SSatish Balay         }
729c4762a1bSJed Brown         case 1: {
730c4762a1bSJed Brown           PetscInt hexSizes_p2[2]  = {5, 4};
731c4762a1bSJed Brown           PetscInt hexPoints_p2[9] = {3, 4, 5, 7, 8, 0, 1, 2, 6};
732c4762a1bSJed Brown 
7339566063dSJacob Faibussowitsch           PetscCall(PetscMalloc2(2, &sizes, 9, &points));
7349566063dSJacob Faibussowitsch           PetscCall(PetscArraycpy(sizes, hexSizes_p2, 2));
7359371c9d4SSatish Balay           PetscCall(PetscArraycpy(points, hexPoints_p2, 9));
7369371c9d4SSatish Balay           break;
7379371c9d4SSatish Balay         }
738d71ae5a4SJacob Faibussowitsch         default:
739d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test number %" PetscInt_FMT " for hexahedral hybrid mesh on 2 procs", user->testNum);
740c4762a1bSJed Brown         }
741c4762a1bSJed Brown       } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Could not find matching test partition");
742c4762a1bSJed Brown     }
7439566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
7449566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetType(part, PETSCPARTITIONERSHELL));
7459566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerShellSetPartition(part, size, sizes, points));
7469566063dSJacob Faibussowitsch     PetscCall(PetscFree2(sizes, points));
747c4762a1bSJed Brown   } else {
748c4762a1bSJed Brown     PetscPartitioner part;
749c4762a1bSJed Brown 
7509566063dSJacob Faibussowitsch     PetscCall(DMPlexGetPartitioner(*dm, &part));
7519566063dSJacob Faibussowitsch     PetscCall(PetscPartitionerSetFromOptions(part));
752c4762a1bSJed Brown   }
753c4762a1bSJed Brown   {
754c4762a1bSJed Brown     DM pdm = NULL;
755c4762a1bSJed Brown 
7569566063dSJacob Faibussowitsch     PetscCall(DMPlexDistribute(*dm, 0, NULL, &pdm));
757c4762a1bSJed Brown     if (pdm) {
7589566063dSJacob Faibussowitsch       PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
7599566063dSJacob Faibussowitsch       PetscCall(DMDestroy(dm));
760c4762a1bSJed Brown       *dm = pdm;
761c4762a1bSJed Brown     }
762c4762a1bSJed Brown   }
7639566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
7649566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view_pre"));
7659566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
766c4762a1bSJed Brown   if (user->uninterpolate || user->reinterpolate) {
767c4762a1bSJed Brown     DM udm = NULL;
768c4762a1bSJed Brown 
7699566063dSJacob Faibussowitsch     PetscCall(DMPlexUninterpolate(*dm, &udm));
7709566063dSJacob Faibussowitsch     PetscCall(DMPlexCopyCoordinates(*dm, udm));
7719566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
772c4762a1bSJed Brown     *dm = udm;
773c4762a1bSJed Brown   }
774c4762a1bSJed Brown   if (user->reinterpolate) {
775c4762a1bSJed Brown     DM idm = NULL;
776c4762a1bSJed Brown 
7779566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
7789566063dSJacob Faibussowitsch     PetscCall(DMPlexCopyCoordinates(*dm, idm));
7799566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
780c4762a1bSJed Brown     *dm = idm;
781c4762a1bSJed Brown   }
7829566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetDefault(*dm, PETSC_FALSE));
7839566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, "Hybrid Mesh"));
7849566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
7859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*dm, "hyb_"));
7869566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(*dm));
7873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
788c4762a1bSJed Brown }
789c4762a1bSJed Brown 
main(int argc,char ** argv)790d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
791d71ae5a4SJacob Faibussowitsch {
792c4762a1bSJed Brown   DM     dm;
793c4762a1bSJed Brown   AppCtx user; /* user-defined work context */
794c4762a1bSJed Brown 
795327415f7SBarry Smith   PetscFunctionBeginUser;
7969566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
7979566063dSJacob Faibussowitsch   PetscCall(ProcessOptions(PETSC_COMM_WORLD, &user));
7989566063dSJacob Faibussowitsch   PetscCall(CreateMesh(PETSC_COMM_WORLD, &user, &dm));
7999566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
8009566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
801b122ec5aSJacob Faibussowitsch   return 0;
802c4762a1bSJed Brown }
803c4762a1bSJed Brown 
804c4762a1bSJed Brown /*TEST
805c4762a1bSJed Brown 
806c4762a1bSJed Brown   # 1D Simplex 29-31
807c4762a1bSJed Brown   testset:
80854fcfd0cSMatthew G. Knepley     args: -dim 1 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
809*3886731fSPierre Jolivet     output_file: output/empty.out
810c4762a1bSJed Brown     test:
811c4762a1bSJed Brown       suffix: 29
812c4762a1bSJed Brown     test:
813c4762a1bSJed Brown       suffix: 30
814c4762a1bSJed Brown       args: -dm_refine 1
815c4762a1bSJed Brown     test:
816c4762a1bSJed Brown       suffix: 31
817c4762a1bSJed Brown       args: -dm_refine 5
818c4762a1bSJed Brown 
819c4762a1bSJed Brown   # 2D Simplex 0-3
820c4762a1bSJed Brown   testset:
82154fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
822*3886731fSPierre Jolivet     output_file: output/empty.out
823c4762a1bSJed Brown     test:
824c4762a1bSJed Brown       suffix: 0
825c4762a1bSJed Brown     test:
826c4762a1bSJed Brown       suffix: 1
8277f9d8d6cSVaclav Hapla       args: -dm_refine 1
828c4762a1bSJed Brown     test:
829c4762a1bSJed Brown       suffix: 2
830c4762a1bSJed Brown       nsize: 2
831c4762a1bSJed Brown     test:
832c4762a1bSJed Brown       suffix: 3
833c4762a1bSJed Brown       nsize: 2
8347f9d8d6cSVaclav Hapla       args: -dm_refine 1
835c4762a1bSJed Brown     test:
836c4762a1bSJed Brown       suffix: 32
837c4762a1bSJed Brown       args: -dm_refine 1 -uninterpolate
838c4762a1bSJed Brown     test:
839c4762a1bSJed Brown       suffix: 33
840c4762a1bSJed Brown       nsize: 2
841c4762a1bSJed Brown       args: -dm_refine 1 -uninterpolate
842c4762a1bSJed Brown     test:
843c4762a1bSJed Brown       suffix: 34
844c4762a1bSJed Brown       nsize: 2
845c4762a1bSJed Brown       args: -dm_refine 3 -uninterpolate
846c4762a1bSJed Brown 
847c4762a1bSJed Brown   # 2D Hybrid Simplex 4-7
848c4762a1bSJed Brown   testset:
84954fcfd0cSMatthew G. Knepley     args: -dim 2 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
850*3886731fSPierre Jolivet     output_file: output/empty.out
851c4762a1bSJed Brown     test:
852c4762a1bSJed Brown       suffix: 4
853c4762a1bSJed Brown     test:
854c4762a1bSJed Brown       suffix: 5
855c4762a1bSJed Brown       args: -dm_refine 1
856c4762a1bSJed Brown     test:
857c4762a1bSJed Brown       suffix: 6
858c4762a1bSJed Brown       nsize: 2
859c4762a1bSJed Brown     test:
860c4762a1bSJed Brown       suffix: 7
861c4762a1bSJed Brown       nsize: 2
862c4762a1bSJed Brown       args: -dm_refine 1
863c4762a1bSJed Brown     test:
864c4762a1bSJed Brown       suffix: 24
865c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
866c4762a1bSJed Brown 
867c4762a1bSJed Brown   # 2D Quad 12-13
868c4762a1bSJed Brown   testset:
86954fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
870*3886731fSPierre Jolivet     output_file: output/empty.out
871c4762a1bSJed Brown     test:
872c4762a1bSJed Brown       suffix: 12
873c4762a1bSJed Brown       args: -dm_refine 1
874c4762a1bSJed Brown     test:
875c4762a1bSJed Brown       suffix: 13
876c4762a1bSJed Brown       nsize: 2
877c4762a1bSJed Brown       args: -dm_refine 1
878c4762a1bSJed Brown 
879c4762a1bSJed Brown   # 2D Hybrid Quad 27-28
880c4762a1bSJed Brown   testset:
88154fcfd0cSMatthew G. Knepley     args: -dim 2 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
882*3886731fSPierre Jolivet     output_file: output/empty.out
883c4762a1bSJed Brown     test:
884c4762a1bSJed Brown       suffix: 27
885c4762a1bSJed Brown       args: -dm_refine 1
886c4762a1bSJed Brown     test:
887c4762a1bSJed Brown       suffix: 28
888c4762a1bSJed Brown       nsize: 2
889c4762a1bSJed Brown       args: -dm_refine 1
890c4762a1bSJed Brown 
891c4762a1bSJed Brown   # 3D Simplex 8-11
892c4762a1bSJed Brown   testset:
89354fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
894*3886731fSPierre Jolivet     output_file: output/empty.out
895c4762a1bSJed Brown     test:
896c4762a1bSJed Brown       suffix: 8
897c4762a1bSJed Brown       args: -dm_refine 1
898c4762a1bSJed Brown     test:
899c4762a1bSJed Brown       suffix: 9
900c4762a1bSJed Brown       nsize: 2
901c4762a1bSJed Brown       args: -dm_refine 1
902c4762a1bSJed Brown     test:
903c4762a1bSJed Brown       suffix: 10
904c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
905c4762a1bSJed Brown     test:
906c4762a1bSJed Brown       suffix: 11
907c4762a1bSJed Brown       nsize: 2
908c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
909c4762a1bSJed Brown     test:
910c4762a1bSJed Brown       suffix: 25
911c4762a1bSJed Brown       args: -test_num 2 -dm_refine 2
912c4762a1bSJed Brown 
913c4762a1bSJed Brown   # 3D Hybrid Simplex 16-19
914c4762a1bSJed Brown   testset:
91554fcfd0cSMatthew G. Knepley     args: -dim 3 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
916*3886731fSPierre Jolivet     output_file: output/empty.out
917c4762a1bSJed Brown     test:
918c4762a1bSJed Brown       suffix: 16
919c4762a1bSJed Brown       args: -dm_refine 1
920c4762a1bSJed Brown     test:
921c4762a1bSJed Brown       suffix: 17
922c4762a1bSJed Brown       nsize: 2
923c4762a1bSJed Brown       args: -dm_refine 1
924c4762a1bSJed Brown     test:
925c4762a1bSJed Brown       suffix: 18
926c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
927c4762a1bSJed Brown     test:
928c4762a1bSJed Brown       suffix: 19
929c4762a1bSJed Brown       nsize: 2
930c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
931c4762a1bSJed Brown 
932c4762a1bSJed Brown   # 3D Hex 14-15
933c4762a1bSJed Brown   testset:
93454fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_simplex 0 -cell_hybrid 0 -hyb_dm_plex_check_all -dm_plex_check_all
935*3886731fSPierre Jolivet     output_file: output/empty.out
936c4762a1bSJed Brown     test:
937c4762a1bSJed Brown       suffix: 14
938c4762a1bSJed Brown       args: -dm_refine 1
939c4762a1bSJed Brown     test:
940c4762a1bSJed Brown       suffix: 15
941c4762a1bSJed Brown       nsize: 2
942c4762a1bSJed Brown      args: -dm_refine 1
943c4762a1bSJed Brown     test:
944c4762a1bSJed Brown       suffix: 26
945c4762a1bSJed Brown       args: -test_num 1 -dm_refine 2
946c4762a1bSJed Brown 
947c4762a1bSJed Brown   # 3D Hybrid Hex 20-23
948c4762a1bSJed Brown   testset:
94954fcfd0cSMatthew G. Knepley     args: -dim 3 -cell_simplex 0 -hyb_dm_plex_check_all -in_dm_plex_check_all -dm_plex_check_all
950*3886731fSPierre Jolivet     output_file: output/empty.out
951c4762a1bSJed Brown     test:
952c4762a1bSJed Brown       suffix: 20
953c4762a1bSJed Brown       args: -dm_refine 1
954c4762a1bSJed Brown     test:
955c4762a1bSJed Brown       suffix: 21
956c4762a1bSJed Brown       nsize: 2
957c4762a1bSJed Brown       args: -dm_refine 1
958c4762a1bSJed Brown     test:
959c4762a1bSJed Brown       suffix: 22
960c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
961c4762a1bSJed Brown     test:
962c4762a1bSJed Brown       suffix: 23
963c4762a1bSJed Brown       nsize: 2
964c4762a1bSJed Brown       args: -test_num 1 -dm_refine 1
965c4762a1bSJed Brown 
966c4762a1bSJed Brown   # Hybrid interpolation
96754fcfd0cSMatthew G. Knepley   #   TODO Setup new tests (like -reinterpolate) that interpolate hybrid cells
968c4762a1bSJed Brown   testset:
969c4762a1bSJed Brown     nsize: 2
97054fcfd0cSMatthew 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
971c4762a1bSJed Brown     test:
972c4762a1bSJed Brown       suffix: hybint_2d_0
973c4762a1bSJed Brown       args: -dim 2 -dm_refine 2
974c4762a1bSJed Brown     test:
975c4762a1bSJed Brown       suffix: hybint_2d_1
976c4762a1bSJed Brown       args: -dim 2 -dm_refine 2 -test_num 1
977c4762a1bSJed Brown     test:
978c4762a1bSJed Brown       suffix: hybint_3d_0
979c4762a1bSJed Brown       args: -dim 3 -dm_refine 1
980c4762a1bSJed Brown     test:
981c4762a1bSJed Brown       suffix: hybint_3d_1
982c4762a1bSJed Brown       args: -dim 3 -dm_refine 1 -test_num 1
983c4762a1bSJed Brown 
984c4762a1bSJed Brown TEST*/
985