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