xref: /petsc/src/dm/impls/plex/tests/ex34.c (revision 2c71b3e237ead271e4f3aa1505f92bf476e3413d)
1c4762a1bSJed Brown static char help[] = "Tests interpolation and output of hybrid meshes\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdmplex.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown typedef struct {
6c4762a1bSJed Brown   char      filename[PETSC_MAX_PATH_LEN]; /* Import mesh from file */
7c4762a1bSJed Brown   PetscBool interpolate;                  /* Interpolate the mesh */
8c4762a1bSJed Brown   PetscInt  meshNum;                      /* Which mesh we should construct */
9c4762a1bSJed Brown } AppCtx;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   PetscErrorCode ierr;
14c4762a1bSJed Brown 
15c4762a1bSJed Brown   PetscFunctionBeginUser;
16c4762a1bSJed Brown   options->filename[0] = '\0';
17c4762a1bSJed Brown   options->interpolate = PETSC_FALSE;
18c4762a1bSJed Brown   options->meshNum     = 0;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown   ierr = PetscOptionsBegin(comm, "", "Hybrid Output Test Options", "DMPLEX");CHKERRQ(ierr);
21589a23caSBarry Smith   ierr = PetscOptionsString("-filename", "The mesh file", "ex8.c", options->filename, options->filename, sizeof(options->filename), NULL);CHKERRQ(ierr);
22c4762a1bSJed Brown   ierr = PetscOptionsBool("-interpolate", "Interpolate the mesh", "ex8.c", options->interpolate, &options->interpolate, NULL);CHKERRQ(ierr);
23c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-mesh_num", "The mesh we should construct", "ex8.c", options->meshNum, &options->meshNum, NULL,0);CHKERRQ(ierr);
241e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   PetscFunctionReturn(0);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown static PetscErrorCode CreateHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   PetscInt       dim;
32c4762a1bSJed Brown   PetscErrorCode ierr;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   PetscFunctionBegin;
35c4762a1bSJed Brown   dim  = 3;
36c4762a1bSJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
37c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Simple Hybrid Mesh");CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
40c4762a1bSJed Brown   {
41c4762a1bSJed Brown     /* Simple mesh with 2 tets and 1 wedge */
42c4762a1bSJed Brown     PetscInt    numPoints[2]         = {8, 3};
43c4762a1bSJed Brown     PetscInt    coneSize[11]         = {4, 4, 6,  0, 0, 0, 0, 0, 0, 0, 0};
44c4762a1bSJed Brown     PetscInt    cones[14]            = {4, 5, 6, 3,  7, 9, 8, 10,  4, 5, 6, 7, 8, 9};
45c4762a1bSJed Brown     PetscInt    coneOrientations[14] = {0, 0, 0, 0,  0, 0, 0,  0,  0, 0, 0, 0, 0, 0};
46c4762a1bSJed Brown     PetscScalar vertexCoords[48]     = {-1.0, 1.0, 0.0,
47c4762a1bSJed Brown                                          0.0, 0.0, 0.0,  0.0, 1.0, -1.0,  0.0, 1.0, 1.0,
48c4762a1bSJed Brown                                          1.0, 0.0, 0.0,  1.0, 1.0, -1.0,  1.0, 1.0, 1.0,
49c4762a1bSJed Brown                                          2.0, 1.0, 0.0};
50c4762a1bSJed Brown 
51c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
52c4762a1bSJed Brown     if (interpolate) {
53c4762a1bSJed Brown       DM idm;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
56c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
57c4762a1bSJed Brown       *dm  = idm;
58c4762a1bSJed Brown     }
59c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
64b5a892a1SMatthew G. Knepley /*
65b5a892a1SMatthew G. Knepley    This is not a valid mesh. We need to either change to tensor quad prisms or regular triangular prisms.
66b5a892a1SMatthew G. Knepley 
67b5a892a1SMatthew G. Knepley            10-------16--------20
68b5a892a1SMatthew G. Knepley            /|        |
69b5a892a1SMatthew G. Knepley           / |        |
70b5a892a1SMatthew G. Knepley          /  |        |
71b5a892a1SMatthew G. Knepley         9---|---15   |
72b5a892a1SMatthew G. Knepley        /|   7    |  13--------18
73b5a892a1SMatthew G. Knepley       / |  /     |  /    ____/
74b5a892a1SMatthew G. Knepley      /  | /      | /____/
75b5a892a1SMatthew G. Knepley     8   |/  14---|//---19
76b5a892a1SMatthew G. Knepley     |   6    |  12
77b5a892a1SMatthew G. Knepley     |  /     |  / \
78b5a892a1SMatthew G. Knepley     | /      | /   \__
79b5a892a1SMatthew G. Knepley     |/       |/       \
80b5a892a1SMatthew G. Knepley     5--------11--------17
81b5a892a1SMatthew G. Knepley */
82c4762a1bSJed Brown static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm)
83c4762a1bSJed Brown {
84c4762a1bSJed Brown   PetscInt       dim;
85c4762a1bSJed Brown   PetscErrorCode ierr;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionBegin;
88c4762a1bSJed Brown   dim  = 3;
89c4762a1bSJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Reverse Hybrid Mesh");CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
93c4762a1bSJed Brown   {
94c4762a1bSJed Brown     /* Simple mesh with 2 hexes and 3 wedges */
95c4762a1bSJed Brown     PetscInt    numPoints[2]         = {16, 5};
96c4762a1bSJed Brown     PetscInt    coneSize[21]         = {8, 8, 6, 6, 6,
97c4762a1bSJed Brown                                         0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
98c4762a1bSJed Brown     PetscInt    cones[34]            = { 5,  6, 12, 11, 8, 14, 15, 9,
99c4762a1bSJed Brown                                          6,  7, 13, 12, 9, 15, 16, 10,
100c4762a1bSJed Brown                                         11, 17, 12, 14, 19, 15,
101c4762a1bSJed Brown                                         12, 18, 13, 15, 20, 16,
102c4762a1bSJed Brown                                         12, 17, 18, 15, 19, 20};
103c4762a1bSJed Brown     PetscInt    coneOrientations[34] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0, 0, 0,
104c4762a1bSJed Brown                                         0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0};
105c4762a1bSJed Brown     PetscScalar vertexCoords[48]     = {-1.0, -1.0, 0.0,  -1.0,  0.0, 0.0,  -1.0, 1.0, 0.0,
106c4762a1bSJed Brown                                         -1.0, -1.0, 1.0,  -1.0,  0.0, 1.0,  -1.0, 1.0, 1.0,
107c4762a1bSJed Brown                                          0.0, -1.0, 0.0,   0.0,  0.0, 0.0,   0.0, 1.0, 0.0,
108c4762a1bSJed Brown                                          0.0, -1.0, 1.0,   0.0,  0.0, 1.0,   0.0, 1.0, 1.0,
109c4762a1bSJed Brown                                          1.0, -1.0, 0.0,                     1.0, 1.0, 0.0,
110c4762a1bSJed Brown                                          1.0, -1.0, 1.0,                     1.0, 1.0, 1.0};
111c4762a1bSJed Brown 
112c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
113c4762a1bSJed Brown     if (interpolate) {
114c4762a1bSJed Brown       DM idm;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
117c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
118c4762a1bSJed Brown       *dm  = idm;
119c4762a1bSJed Brown     }
120c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
121c4762a1bSJed Brown   }
122c4762a1bSJed Brown   PetscFunctionReturn(0);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown static PetscErrorCode OrderHybridMesh(DM *dm)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   DM             pdm;
128c4762a1bSJed Brown   IS             perm;
129c4762a1bSJed Brown   PetscInt      *ind;
130c4762a1bSJed Brown   PetscInt       dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2];
131c4762a1bSJed Brown   PetscErrorCode ierr;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBegin;
134c4762a1bSJed Brown   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
135*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(dim != 3,PetscObjectComm((PetscObject) *dm), PETSC_ERR_SUP, "No support for dimension %D", dim);
136c4762a1bSJed Brown   ierr = DMPlexGetChart(*dm, &pStart, &pEnd);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = PetscMalloc1(pEnd-pStart, &ind);CHKERRQ(ierr);
138c4762a1bSJed Brown   for (p = 0; p < pEnd-pStart; ++p) ind[p] = p;
139c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
140c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
141c4762a1bSJed Brown     PetscInt coneSize;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown     ierr = DMPlexGetConeSize(*dm, c, &coneSize);CHKERRQ(ierr);
144c4762a1bSJed Brown     if (coneSize == 6) ++Nhyb;
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown   off[0] = 0;
147c4762a1bSJed Brown   off[1] = cEnd - Nhyb;
148c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
149c4762a1bSJed Brown     PetscInt coneSize;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown     ierr = DMPlexGetConeSize(*dm, c, &coneSize);CHKERRQ(ierr);
152c4762a1bSJed Brown     if (coneSize == 6) ind[c] = off[1]++;
153c4762a1bSJed Brown     else               ind[c] = off[0]++;
154c4762a1bSJed Brown   }
155*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(off[0] != cEnd - Nhyb,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of normal cells %D should be %D", off[0], cEnd - Nhyb);
156*2c71b3e2SJacob Faibussowitsch   PetscCheckFalse(off[1] != cEnd,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of hybrid cells %D should be %D", off[1] - off[0], Nhyb);
157c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, ind, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
158c4762a1bSJed Brown   ierr = DMPlexPermute(*dm, perm, &pdm);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
160c4762a1bSJed Brown   ierr = DMDestroy(dm);CHKERRQ(ierr);
161c4762a1bSJed Brown   *dm  = pdm;
162c4762a1bSJed Brown   PetscFunctionReturn(0);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
165c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
166c4762a1bSJed Brown {
167c4762a1bSJed Brown   const char    *filename    = user->filename;
168c4762a1bSJed Brown   PetscBool      interpolate = user->interpolate;
169c4762a1bSJed Brown   PetscInt       meshNum     = user->meshNum;
170c4762a1bSJed Brown   size_t         len;
171c4762a1bSJed Brown   PetscErrorCode ierr;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   PetscFunctionBegin;
174c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
175c4762a1bSJed Brown   if (len) {
176cd7e8a5eSksagiyam     ierr = DMPlexCreateFromFile(comm, filename, "ex34_plex", PETSC_FALSE, dm);CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = OrderHybridMesh(dm);CHKERRQ(ierr);
178c4762a1bSJed Brown     if (interpolate) {
179c4762a1bSJed Brown       DM idm;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
182c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
183c4762a1bSJed Brown       *dm  = idm;
184c4762a1bSJed Brown     }
185c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) *dm, "Input Mesh");CHKERRQ(ierr);
186c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
187c4762a1bSJed Brown   } else {
188c4762a1bSJed Brown     switch (meshNum) {
189c4762a1bSJed Brown     case 0:
190c4762a1bSJed Brown       ierr = CreateHybridMesh(comm, interpolate, dm);CHKERRQ(ierr);break;
191c4762a1bSJed Brown     case 1:
192c4762a1bSJed Brown       ierr = CreateReverseHybridMesh(comm, interpolate, dm);CHKERRQ(ierr);break;
19398921bdaSJacob Faibussowitsch     default: SETERRQ(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %D", user->meshNum);
194c4762a1bSJed Brown     }
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown int main(int argc, char **argv)
200c4762a1bSJed Brown {
201c4762a1bSJed Brown   DM             dm;
202c4762a1bSJed Brown   AppCtx         user;
203c4762a1bSJed Brown   PetscErrorCode ierr;
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
206c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
207c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
208c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
209c4762a1bSJed Brown   ierr = PetscFinalize();
210c4762a1bSJed Brown   return ierr;
211c4762a1bSJed Brown }
212c4762a1bSJed Brown 
213c4762a1bSJed Brown /*TEST
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   test:
216c4762a1bSJed Brown     suffix: 0
217c4762a1bSJed Brown     args: -interpolate -dm_view ascii::ascii_info_detail
218c4762a1bSJed Brown 
219b5a892a1SMatthew G. Knepley   # Test needs to be reworked
220c4762a1bSJed Brown   test:
221b5a892a1SMatthew G. Knepley     requires: BROKEN
222c4762a1bSJed Brown     suffix: 1
223c4762a1bSJed Brown     args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail
224c4762a1bSJed Brown 
225c4762a1bSJed Brown TEST*/
226