xref: /petsc/src/dm/impls/plex/tests/ex34.c (revision 589a23caa660d2a5f330cc8d1ed213e9cfaf51a7)
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);
21*589a23caSBarry 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);
24c4762a1bSJed Brown   ierr = PetscOptionsEnd();
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 
64c4762a1bSJed Brown /* TODO Need to extend interpolation to work when regular cells join to hybrid faces, rather than end faces */
65c4762a1bSJed Brown static PetscErrorCode CreateReverseHybridMesh(MPI_Comm comm, PetscBool interpolate, DM *dm)
66c4762a1bSJed Brown {
67c4762a1bSJed Brown   PetscInt       dim;
68c4762a1bSJed Brown   PetscErrorCode ierr;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBegin;
71c4762a1bSJed Brown   dim  = 3;
72c4762a1bSJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
73c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, "Reverse Hybrid Mesh");CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
76c4762a1bSJed Brown   {
77c4762a1bSJed Brown     /* Simple mesh with 2 hexes and 3 wedges */
78c4762a1bSJed Brown     PetscInt    numPoints[2]         = {16, 5};
79c4762a1bSJed Brown     PetscInt    coneSize[21]         = {8, 8, 6, 6, 6,
80c4762a1bSJed Brown                                         0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0};
81c4762a1bSJed Brown     PetscInt    cones[34]            = { 5,  6, 12, 11, 8, 14, 15, 9,
82c4762a1bSJed Brown                                          6,  7, 13, 12, 9, 15, 16, 10,
83c4762a1bSJed Brown                                         11, 17, 12, 14, 19, 15,
84c4762a1bSJed Brown                                         12, 18, 13, 15, 20, 16,
85c4762a1bSJed Brown                                         12, 17, 18, 15, 19, 20};
86c4762a1bSJed Brown     PetscInt    coneOrientations[34] = {0, 0, 0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0, 0, 0,
87c4762a1bSJed Brown                                         0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0,  0, 0, 0, 0, 0, 0};
88c4762a1bSJed Brown     PetscScalar vertexCoords[48]     = {-1.0, -1.0, 0.0,  -1.0,  0.0, 0.0,  -1.0, 1.0, 0.0,
89c4762a1bSJed Brown                                         -1.0, -1.0, 1.0,  -1.0,  0.0, 1.0,  -1.0, 1.0, 1.0,
90c4762a1bSJed Brown                                          0.0, -1.0, 0.0,   0.0,  0.0, 0.0,   0.0, 1.0, 0.0,
91c4762a1bSJed Brown                                          0.0, -1.0, 1.0,   0.0,  0.0, 1.0,   0.0, 1.0, 1.0,
92c4762a1bSJed Brown                                          1.0, -1.0, 0.0,                     1.0, 1.0, 0.0,
93c4762a1bSJed Brown                                          1.0, -1.0, 1.0,                     1.0, 1.0, 1.0};
94c4762a1bSJed Brown 
95c4762a1bSJed Brown     ierr = DMPlexCreateFromDAG(*dm, 1, numPoints, coneSize, cones, coneOrientations, vertexCoords);CHKERRQ(ierr);
96c4762a1bSJed Brown     if (interpolate) {
97c4762a1bSJed Brown       DM idm;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
100c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
101c4762a1bSJed Brown       *dm  = idm;
102c4762a1bSJed Brown     }
103c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown   PetscFunctionReturn(0);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown static PetscErrorCode OrderHybridMesh(DM *dm)
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   DM             pdm;
111c4762a1bSJed Brown   IS             perm;
112c4762a1bSJed Brown   PetscInt      *ind;
113c4762a1bSJed Brown   PetscInt       dim, pStart, pEnd, p, cStart, cEnd, c, Nhyb = 0, off[2];
114c4762a1bSJed Brown   PetscErrorCode ierr;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBegin;
117c4762a1bSJed Brown   ierr = DMGetDimension(*dm, &dim);CHKERRQ(ierr);
118c4762a1bSJed Brown   if (dim != 3) SETERRQ1(PetscObjectComm((PetscObject) *dm), PETSC_ERR_SUP, "No support for dimension %D", dim);
119c4762a1bSJed Brown   ierr = DMPlexGetChart(*dm, &pStart, &pEnd);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = PetscMalloc1(pEnd-pStart, &ind);CHKERRQ(ierr);
121c4762a1bSJed Brown   for (p = 0; p < pEnd-pStart; ++p) ind[p] = p;
122c4762a1bSJed Brown   ierr = DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
123c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
124c4762a1bSJed Brown     PetscInt coneSize;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown     ierr = DMPlexGetConeSize(*dm, c, &coneSize);CHKERRQ(ierr);
127c4762a1bSJed Brown     if (coneSize == 6) ++Nhyb;
128c4762a1bSJed Brown   }
129c4762a1bSJed Brown   off[0] = 0;
130c4762a1bSJed Brown   off[1] = cEnd - Nhyb;
131c4762a1bSJed Brown   for (c = cStart; c < cEnd; ++c) {
132c4762a1bSJed Brown     PetscInt coneSize;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown     ierr = DMPlexGetConeSize(*dm, c, &coneSize);CHKERRQ(ierr);
135c4762a1bSJed Brown     if (coneSize == 6) ind[c] = off[1]++;
136c4762a1bSJed Brown     else               ind[c] = off[0]++;
137c4762a1bSJed Brown   }
138c4762a1bSJed Brown   if (off[0] != cEnd - Nhyb) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of normal cells %D should be %D", off[0], cEnd - Nhyb);
139c4762a1bSJed Brown   if (off[1] != cEnd)        SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of hybrid cells %D should be %D", off[1] - off[0], Nhyb);
140c4762a1bSJed Brown   ierr = ISCreateGeneral(PETSC_COMM_SELF, pEnd-pStart, ind, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr = DMPlexPermute(*dm, perm, &pdm);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr = DMDestroy(dm);CHKERRQ(ierr);
144c4762a1bSJed Brown   *dm  = pdm;
145c4762a1bSJed Brown   PetscFunctionReturn(0);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *user, DM *dm)
149c4762a1bSJed Brown {
150c4762a1bSJed Brown   const char    *filename    = user->filename;
151c4762a1bSJed Brown   PetscBool      interpolate = user->interpolate;
152c4762a1bSJed Brown   PetscInt       meshNum     = user->meshNum;
153c4762a1bSJed Brown   size_t         len;
154c4762a1bSJed Brown   PetscErrorCode ierr;
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   PetscFunctionBegin;
157c4762a1bSJed Brown   ierr = PetscStrlen(filename, &len);CHKERRQ(ierr);
158c4762a1bSJed Brown   if (len) {
159c4762a1bSJed Brown     ierr = DMPlexCreateFromFile(comm, filename, PETSC_FALSE, dm);CHKERRQ(ierr);
160c4762a1bSJed Brown     ierr = OrderHybridMesh(dm);CHKERRQ(ierr);
161c4762a1bSJed Brown     if (interpolate) {
162c4762a1bSJed Brown       DM idm;
163c4762a1bSJed Brown 
164c4762a1bSJed Brown       ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
165c4762a1bSJed Brown       ierr = DMDestroy(dm);CHKERRQ(ierr);
166c4762a1bSJed Brown       *dm  = idm;
167c4762a1bSJed Brown     }
168c4762a1bSJed Brown     ierr = PetscObjectSetName((PetscObject) *dm, "Input Mesh");CHKERRQ(ierr);
169c4762a1bSJed Brown     ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
170c4762a1bSJed Brown   } else {
171c4762a1bSJed Brown     switch (meshNum) {
172c4762a1bSJed Brown     case 0:
173c4762a1bSJed Brown       ierr = CreateHybridMesh(comm, interpolate, dm);CHKERRQ(ierr);break;
174c4762a1bSJed Brown     case 1:
175c4762a1bSJed Brown       ierr = CreateReverseHybridMesh(comm, interpolate, dm);CHKERRQ(ierr);break;
176c4762a1bSJed Brown     default: SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "Unknown mesh number %D", user->meshNum);
177c4762a1bSJed Brown     }
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   PetscFunctionReturn(0);
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown int main(int argc, char **argv)
183c4762a1bSJed Brown {
184c4762a1bSJed Brown   DM             dm;
185c4762a1bSJed Brown   AppCtx         user;
186c4762a1bSJed Brown   PetscErrorCode ierr;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
189c4762a1bSJed Brown   ierr = ProcessOptions(PETSC_COMM_WORLD, &user);CHKERRQ(ierr);
190c4762a1bSJed Brown   ierr = CreateMesh(PETSC_COMM_WORLD, &user, &dm);CHKERRQ(ierr);
191c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
192c4762a1bSJed Brown   ierr = PetscFinalize();
193c4762a1bSJed Brown   return ierr;
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /*TEST
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   test:
199c4762a1bSJed Brown     suffix: 0
200c4762a1bSJed Brown     args: -interpolate -dm_view ascii::ascii_info_detail
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   test:
203c4762a1bSJed Brown     suffix: 1
204c4762a1bSJed Brown     args: -mesh_num 1 -interpolate -dm_view ascii::ascii_info_detail
205c4762a1bSJed Brown 
206c4762a1bSJed Brown TEST*/
207