xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision df6f6c19b095f1abc1ac27fcbb9c2092a4a2a65e)
1 #define PETSCDM_DLL
2 #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMPlexCreateGmshFromFile"
6 /*@C
7   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
8 
9 + comm        - The MPI communicator
10 . filename    - Name of the Gmsh file
11 - interpolate - Create faces and edges in the mesh
12 
13   Output Parameter:
14 . dm  - The DM object representing the mesh
15 
16   Level: beginner
17 
18 .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
19 @*/
20 PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
21 {
22   PetscViewer     viewer, vheader;
23   PetscMPIInt     rank;
24   PetscViewerType vtype;
25   char            line[PETSC_MAX_PATH_LEN];
26   int             snum;
27   PetscBool       match;
28   int             fileType;
29   PetscErrorCode  ierr;
30 
31   PetscFunctionBegin;
32   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
33   /* Determine Gmsh file type (ASCII or binary) from file header */
34   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
35   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
36   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
37   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
38   if (!rank) {
39     /* Read only the first two lines of the Gmsh file */
40     ierr = PetscViewerRead(vheader, line, 1, PETSC_STRING);CHKERRQ(ierr);
41     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
42     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
43     ierr = PetscViewerRead(vheader, line, 2, PETSC_STRING);CHKERRQ(ierr);
44     snum = sscanf(line, "2.2 %d", &fileType);
45     if (snum != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
46   }
47   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
48   /* Create appropriate viewer and build plex */
49   if (fileType == 0) vtype = PETSCVIEWERASCII;
50   else vtype = PETSCVIEWERBINARY;
51   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
52   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
53   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
54   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
55   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
56   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
57   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
58   PetscFunctionReturn(0);
59 }
60 
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "DMPlexCreateGmsh"
64 /*@
65   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
66 
67   Collective on comm
68 
69   Input Parameters:
70 + comm  - The MPI communicator
71 . viewer - The Viewer associated with a Gmsh file
72 - interpolate - Create faces and edges in the mesh
73 
74   Output Parameter:
75 . dm  - The DM object representing the mesh
76 
77   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
78 
79   Level: beginner
80 
81 .keywords: mesh,Gmsh
82 .seealso: DMPLEX, DMCreate()
83 @*/
84 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
85 {
86   PetscViewerType vtype;
87   GmshElement   *gmsh_elem;
88   PetscSection   coordSection;
89   Vec            coordinates;
90   PetscScalar   *coords, *coordsIn = NULL;
91   PetscInt       dim = 0, coordSize, c, v, d, cell;
92   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
93   PetscMPIInt    num_proc, rank;
94   char           line[PETSC_MAX_PATH_LEN];
95   PetscBool      match, binary, bswap = PETSC_FALSE;
96   PetscErrorCode ierr;
97 
98   PetscFunctionBegin;
99   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
100   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
101   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
102   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
103   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
104   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
105   if (!rank || binary) {
106     PetscBool match;
107     int       fileType, dataSize;
108 
109     /* Read header */
110     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
111     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
112     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
113     ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr);
114     snum = sscanf(line, "2.2 %d %d", &fileType, &dataSize);
115     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
116     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
117     if (binary) {
118       int checkInt;
119       ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_ENUM);CHKERRQ(ierr);
120       if (checkInt != 1) {
121         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
122         if (checkInt == 1) bswap = PETSC_TRUE;
123         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
124       }
125     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
126     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
127     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
128     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
129     /* Read vertices */
130     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
131     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
132     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
133     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
134     snum = sscanf(line, "%d", &numVertices);
135     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
136     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
137     for (v = 0; v < numVertices; ++v) {
138       int i;
139       ierr = PetscViewerRead(viewer, &i, 1, PETSC_ENUM);CHKERRQ(ierr);
140       if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr);
141       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr);
142       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
143       if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
144     }
145     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
146     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
147     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
148     /* Read cells */
149     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
150     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
151     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
152     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
153     snum = sscanf(line, "%d", &numCells);
154     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
155   }
156 
157   if (!rank || binary) {
158     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
159        file contents multiple times to figure out the true number of cells and facets
160        in the given mesh. To make this more efficient we read the file contents only
161        once and store them in memory, while determining the true number of cells. */
162     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
163     for (trueNumCells=0, c = 0; c < numCells; ++c) {
164       ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr);
165       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
166       if (gmsh_elem[c].dim == dim) trueNumCells++;
167     }
168     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
169     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
170     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
171   }
172   /* For binary we read on all ranks, but only build the plex on rank 0 */
173   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
174   /* Allocate the cell-vertex mesh */
175   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
176   if (!rank) {
177     for (cell = 0, c = 0; c < numCells; ++c) {
178       if (gmsh_elem[c].dim == dim) {
179         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
180         cell++;
181       }
182     }
183   }
184   ierr = DMSetUp(*dm);CHKERRQ(ierr);
185   /* Add cell-vertex connections */
186   if (!rank) {
187     PetscInt pcone[8], corner;
188     for (cell = 0, c = 0; c < numCells; ++c) {
189       if (gmsh_elem[c].dim == dim) {
190         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
191           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
192         }
193         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
194         cell++;
195       }
196     }
197   }
198   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
199   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
200   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
201   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
202   if (interpolate) {
203     DM idm = NULL;
204 
205     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
206     ierr = DMDestroy(dm);CHKERRQ(ierr);
207     *dm  = idm;
208   }
209 
210   if (!rank) {
211     /* Apply boundary IDs by finding the relevant facets with vertex joins */
212     PetscInt pcone[8], corner, vStart, vEnd;
213 
214     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
215     for (c = 0; c < numCells; ++c) {
216       if (gmsh_elem[c].dim == dim-1) {
217         PetscInt joinSize;
218         const PetscInt *join;
219         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
220           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
221         }
222         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
223         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
224         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
225         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
226       }
227     }
228   }
229 
230   /* Read coordinates */
231   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
232   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
233   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
234   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
235   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
236     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
237     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
238   }
239   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
240   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
241   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
242   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
243   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
244   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
245   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
246   if (!rank) {
247     for (v = 0; v < numVertices; ++v) {
248       for (d = 0; d < dim; ++d) {
249         coords[v*dim+d] = coordsIn[v*3+d];
250       }
251     }
252   }
253   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
254   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
255   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
256   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
257   /* Clean up intermediate storage */
258   if (!rank || binary) {
259     for (c = 0; c < numCells; ++c) {
260       ierr = PetscFree(gmsh_elem[c].nodes);
261       ierr = PetscFree(gmsh_elem[c].tags);
262     }
263     ierr = PetscFree(gmsh_elem);
264   }
265   PetscFunctionReturn(0);
266 }
267 
268 #undef __FUNCT__
269 #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
270 PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele)
271 {
272   int            cellType, numElem;
273   PetscErrorCode ierr;
274 
275   PetscFunctionBegin;
276   if (binary) {
277     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
278     if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_ENUM, 1);CHKERRQ(ierr);
279     ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_ENUM);CHKERRQ(ierr);
280     if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_ENUM, 1);CHKERRQ(ierr);
281     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
282     if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_ENUM, 1);CHKERRQ(ierr);
283     ierr = PetscViewerRead(viewer,  &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
284     if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_ENUM, 1);CHKERRQ(ierr);
285   } else {
286     ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
287     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
288     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
289   }
290   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
291   ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_ENUM);CHKERRQ(ierr);
292   if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_ENUM, ele->numTags);CHKERRQ(ierr);
293   switch (cellType) {
294   case 1: /* 2-node line */
295     ele->dim = 1;
296     ele->numNodes = 2;
297     break;
298   case 2: /* 3-node triangle */
299     ele->dim = 2;
300     ele->numNodes = 3;
301     break;
302   case 3: /* 4-node quadrangle */
303     ele->dim = 2;
304     ele->numNodes = 4;
305     break;
306   case 4: /* 4-node tetrahedron */
307     ele->dim  = 3;
308     ele->numNodes = 4;
309     break;
310   case 5: /* 8-node hexahedron */
311     ele->dim = 3;
312     ele->numNodes = 8;
313     break;
314   case 15: /* 1-node vertex */
315     ele->dim = 0;
316     ele->numNodes = 1;
317     break;
318   default:
319     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
320   }
321   ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
322   ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_ENUM);CHKERRQ(ierr);
323   if (byteSwap) {ierr = PetscByteSwap(ele->nodes, PETSC_ENUM, ele->numNodes);CHKERRQ(ierr);}
324   PetscFunctionReturn(0);
325 }
326