xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision ecdc758bc49a6338bb732f704b13831ad7419275)
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   FILE          *fd;
29   PetscSection   coordSection;
30   Vec            coordinates;
31   PetscScalar   *coords, *coordsIn = NULL;
32   PetscInt       dim = 0, tdim = 0, coordSize, c, v, d, numCorners, cell;
33   int            numVertices = 0, numCells = 0, trueNumCells = 0, cone[8], tags[2], cellNum, snum;
34   long           fpos = 0;
35   PetscMPIInt    num_proc, rank;
36   char           line[PETSC_MAX_PATH_LEN];
37   PetscBool      match;
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   if (!rank) {
46     PetscBool match;
47     int       fileType, dataSize;
48 
49     ierr = PetscViewerASCIIGetPointer(viewer, &fd);CHKERRQ(ierr);
50     /* Read header */
51     fgets(line, PETSC_MAX_PATH_LEN, fd);
52     ierr = PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
53     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
54     snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2);
55     if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
56     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
57     fgets(line, PETSC_MAX_PATH_LEN, fd);
58     ierr = PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
59     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
60     /* Read vertices */
61     fgets(line, PETSC_MAX_PATH_LEN, fd);
62     ierr = PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
63     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
64     snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1);
65     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
66     for (v = 0; v < numVertices; ++v) {
67       double x, y, z;
68       int    i;
69 
70       snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4);
71       coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z;
72       if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
73     }
74     fgets(line, PETSC_MAX_PATH_LEN, fd);
75     ierr = PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
76     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
77     /* Read cells */
78     fgets(line, PETSC_MAX_PATH_LEN, fd);
79     ierr = PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
80     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
81     snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1);
82   }
83 
84   if (!rank) {
85     fpos = ftell(fd);
86     /* The Gmsh format disguises facets as elements, so we have to run through all "element" entries
87        to get the correct numCells and decide the topological dimension of the mesh */
88     trueNumCells = 0;
89     for (c = 0; c < numCells; ++c) {
90       ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, tags);CHKERRQ(ierr);
91       if (dim > tdim) {
92         tdim = dim;
93         trueNumCells = 0;
94       }
95       if (dim == tdim) trueNumCells++;
96     }
97   }
98   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
99   if (!rank) {
100     ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr);
101     for (cell = 0, c = 0; c < numCells; ++c) {
102       ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, tags);CHKERRQ(ierr);
103       if (dim == tdim) {
104         ierr = DMPlexSetConeSize(*dm, cell, numCorners);CHKERRQ(ierr);
105         cell++;
106       }
107       if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1);
108     }
109   }
110   ierr = DMSetUp(*dm);CHKERRQ(ierr);
111   if (!rank) {
112     PetscInt pcone[8], corner;
113 
114     ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr);
115     for (cell = 0, c = 0; c < numCells; ++c) {
116       ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, tags);CHKERRQ(ierr);
117       if (dim == tdim) {
118         for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + trueNumCells-1;
119         ierr = DMPlexSetCone(*dm, cell, (const PetscInt *) pcone);CHKERRQ(ierr);
120         cell++;
121       }
122       if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1);
123     }
124     fgets(line, PETSC_MAX_PATH_LEN, fd);
125     ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
126     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
127   }
128   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
129   ierr = DMSetDimension(*dm, tdim);CHKERRQ(ierr);
130   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
131   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
132   if (interpolate) {
133     DM idm = NULL;
134 
135     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
136     ierr = DMDestroy(dm);CHKERRQ(ierr);
137     *dm  = idm;
138   }
139 
140   if (!rank) {
141     /* Apply boundary IDs by finding the relevant facets with vertex joins */
142     PetscInt pcone[8], corner, vStart, vEnd;
143 
144     ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr);
145     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
146     for (c = 0; c < numCells; ++c) {
147       ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, tags);CHKERRQ(ierr);
148       if (dim == tdim-1) {
149         PetscInt joinSize;
150         const PetscInt *join;
151         for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + vStart - 1;
152         ierr = DMPlexGetFullJoin(*dm, numCorners, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
153         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", cellNum);
154         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], tags[0]);CHKERRQ(ierr);
155         ierr = DMPlexRestoreJoin(*dm, numCorners, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
156       }
157       if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1);
158     }
159     fgets(line, PETSC_MAX_PATH_LEN, fd);
160     ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
161     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
162   }
163 
164   /* Read coordinates */
165   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
166   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
167   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
168   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
169   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
170     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
171     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
172   }
173   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
174   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
175   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
176   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
177   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
178   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
179   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
180   if (!rank) {
181     for (v = 0; v < numVertices; ++v) {
182       for (d = 0; d < dim; ++d) {
183         coords[v*dim+d] = coordsIn[v*3+d];
184       }
185     }
186   }
187   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
188   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
189   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
190   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
191   PetscFunctionReturn(0);
192 }
193 
194 #undef __FUNCT__
195 #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
196 PetscErrorCode DMPlexCreateGmsh_ReadElement(FILE *fd, PetscInt *dim, int *cellNum, PetscInt *numCorners, int cone[], int tags[])
197 {
198   PetscInt t;
199   int      numTags, snum, cellType;
200 
201   PetscFunctionBegin;
202   snum = fscanf(fd, "%d %d %d", cellNum, &cellType, &numTags);CHKERRQ(snum != 3);
203   for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tags[t]);CHKERRQ(snum != 1);}
204   switch (cellType) {
205   case 1: /* 2-node line */
206     *dim = 1;
207     *numCorners = 2;
208     snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != *numCorners);
209     break;
210   case 2: /* 3-node triangle */
211     *dim = 2;
212     *numCorners = 3;
213     snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != *numCorners);
214     break;
215   case 3: /* 4-node quadrangle */
216     *dim = 2;
217     *numCorners = 4;
218     snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners);
219     break;
220   case 4: /* 4-node tetrahedron */
221     *dim  = 3;
222     *numCorners = 4;
223     snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners);
224     break;
225   case 5: /* 8-node hexahedron */
226     *dim = 3;
227     *numCorners = 8;
228     snum = fscanf(fd, "%d %d %d %d %d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3], &cone[4], &cone[5], &cone[6], &cone[7]);CHKERRQ(snum != *numCorners);
229     break;
230   case 15: /* 1-node vertex */
231     *dim = 0;
232     *numCorners = 1;
233     snum = fscanf(fd, "%d\n", &cone[0]);CHKERRQ(snum != *numCorners);
234     break;
235   default:
236     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
237   }
238   PetscFunctionReturn(0);
239 }
240