xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision a4bb7517d7207e6380c16a59a3e8a6843fbf2279) !
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2331830f3SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3331830f3SMatthew G. Knepley 
4331830f3SMatthew G. Knepley #undef __FUNCT__
5331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
6331830f3SMatthew G. Knepley /*@
7331830f3SMatthew G. Knepley   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file.
8331830f3SMatthew G. Knepley 
9331830f3SMatthew G. Knepley   Collective on comm
10331830f3SMatthew G. Knepley 
11331830f3SMatthew G. Knepley   Input Parameters:
12331830f3SMatthew G. Knepley + comm  - The MPI communicator
13331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
14331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
15331830f3SMatthew G. Knepley 
16331830f3SMatthew G. Knepley   Output Parameter:
17331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
18331830f3SMatthew G. Knepley 
19331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
20331830f3SMatthew G. Knepley 
21331830f3SMatthew G. Knepley   Level: beginner
22331830f3SMatthew G. Knepley 
23331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
24331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
25331830f3SMatthew G. Knepley @*/
26331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
27331830f3SMatthew G. Knepley {
28331830f3SMatthew G. Knepley   FILE          *fd;
29*a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
30331830f3SMatthew G. Knepley   PetscSection   coordSection;
31331830f3SMatthew G. Knepley   Vec            coordinates;
32331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
33*a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
34*a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
35331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
36331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
37331830f3SMatthew G. Knepley   PetscBool      match;
38331830f3SMatthew G. Knepley   PetscErrorCode ierr;
39331830f3SMatthew G. Knepley 
40331830f3SMatthew G. Knepley   PetscFunctionBegin;
41331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
42331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
43331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
44331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
45331830f3SMatthew G. Knepley   if (!rank) {
46331830f3SMatthew G. Knepley     PetscBool match;
47331830f3SMatthew G. Knepley     int       fileType, dataSize;
48331830f3SMatthew G. Knepley 
49331830f3SMatthew G. Knepley     ierr = PetscViewerASCIIGetPointer(viewer, &fd);CHKERRQ(ierr);
50331830f3SMatthew G. Knepley     /* Read header */
51331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
52331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
53331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
54331830f3SMatthew G. Knepley     snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2);
55331830f3SMatthew G. Knepley     if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
56331830f3SMatthew G. Knepley     if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize);
57331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
58331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
59331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
60331830f3SMatthew G. Knepley     /* Read vertices */
61331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
62331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
63331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
64331830f3SMatthew G. Knepley     snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1);
65331830f3SMatthew G. Knepley     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
66331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
67331830f3SMatthew G. Knepley       double x, y, z;
68331830f3SMatthew G. Knepley       int    i;
69331830f3SMatthew G. Knepley 
70331830f3SMatthew G. Knepley       snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4);
71331830f3SMatthew G. Knepley       coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z;
72331830f3SMatthew G. Knepley       if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
73331830f3SMatthew G. Knepley     }
74331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
75331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
76331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
77331830f3SMatthew G. Knepley     /* Read cells */
78331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
79331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
80331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
81331830f3SMatthew G. Knepley     snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1);
82331830f3SMatthew G. Knepley   }
83db397164SMichael Lange 
84331830f3SMatthew G. Knepley   if (!rank) {
85*a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
86*a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
87*a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
88*a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
89*a4bb7517SMichael Lange     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
90*a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
91*a4bb7517SMichael Lange       ierr = DMPlexCreateGmsh_ReadElement(fd, &gmsh_elem[c]);CHKERRQ(ierr);
92*a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
93*a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
94331830f3SMatthew G. Knepley     }
95331830f3SMatthew G. Knepley     fgets(line, PETSC_MAX_PATH_LEN, fd);
96331830f3SMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
97331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
98331830f3SMatthew G. Knepley   }
99*a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
100*a4bb7517SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
101*a4bb7517SMichael Lange   if (!rank) {
102*a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
103*a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
104*a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
105*a4bb7517SMichael Lange         cell++;
106*a4bb7517SMichael Lange       }
107*a4bb7517SMichael Lange     }
108*a4bb7517SMichael Lange   }
109*a4bb7517SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
110*a4bb7517SMichael Lange   /* Add cell-vertex connections */
111*a4bb7517SMichael Lange   if (!rank) {
112*a4bb7517SMichael Lange     PetscInt pcone[8], corner;
113*a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
114*a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
115*a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
116*a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
117*a4bb7517SMichael Lange         }
118*a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
119*a4bb7517SMichael Lange         cell++;
120*a4bb7517SMichael Lange       }
121*a4bb7517SMichael Lange     }
122*a4bb7517SMichael Lange   }
1236228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
124*a4bb7517SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
125331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
126331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
127331830f3SMatthew G. Knepley   if (interpolate) {
128e56d480eSMatthew G. Knepley     DM idm = NULL;
129331830f3SMatthew G. Knepley 
130331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
131331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
132331830f3SMatthew G. Knepley     *dm  = idm;
133331830f3SMatthew G. Knepley   }
13416ad7e67SMichael Lange 
13516ad7e67SMichael Lange   if (!rank) {
13616ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
13716ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
138d1a54cefSMatthew G. Knepley 
13916ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
14016ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
141*a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
14216ad7e67SMichael Lange         PetscInt joinSize;
14316ad7e67SMichael Lange         const PetscInt *join;
144*a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
145*a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
14616ad7e67SMichael Lange         }
147*a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
148*a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
149*a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
150*a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
15116ad7e67SMichael Lange       }
152*a4bb7517SMichael Lange     }
15316ad7e67SMichael Lange   }
15416ad7e67SMichael Lange 
155331830f3SMatthew G. Knepley   /* Read coordinates */
156979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
157331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
158331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
15975b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
16075b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
161331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
162331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
163331830f3SMatthew G. Knepley   }
164331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
165331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
166331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
167331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
168331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
169331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
170331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
171331830f3SMatthew G. Knepley   if (!rank) {
172331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
173331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
174331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
175331830f3SMatthew G. Knepley       }
176331830f3SMatthew G. Knepley     }
177331830f3SMatthew G. Knepley   }
178331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
179331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
180331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
181331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
182*a4bb7517SMichael Lange   /* Clean up intermediate storage */
183*a4bb7517SMichael Lange   if (!rank) {
184*a4bb7517SMichael Lange     for (c = 0; c < numCells; ++c) {
185*a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].nodes);
186*a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].tags);
187*a4bb7517SMichael Lange     }
188*a4bb7517SMichael Lange     ierr = PetscFree(gmsh_elem);
189*a4bb7517SMichael Lange   }
190331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
191331830f3SMatthew G. Knepley }
192db397164SMichael Lange 
193db397164SMichael Lange #undef __FUNCT__
19430412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
195*a4bb7517SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(FILE *fd, GmshElement *ele)
196db397164SMichael Lange {
197*a4bb7517SMichael Lange   int            snum, cellType;
19830412aabSMichael Lange   PetscInt       t;
199*a4bb7517SMichael Lange   PetscErrorCode ierr;
200db397164SMichael Lange 
20130412aabSMichael Lange   PetscFunctionBegin;
202*a4bb7517SMichael Lange   snum = fscanf(fd, "%d %d %d", &(ele->id), &cellType, &(ele->numTags));CHKERRQ(snum != 3);
203*a4bb7517SMichael Lange   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
204*a4bb7517SMichael Lange   for (t=0; t<ele->numTags; t++) {snum = fscanf(fd, "%d", &(ele->tags[t]));CHKERRQ(snum != 1);}
20530412aabSMichael Lange   switch (cellType) {
20630412aabSMichael Lange   case 1: /* 2-node line */
207*a4bb7517SMichael Lange     ele->dim = 1;
208*a4bb7517SMichael Lange     ele->numNodes = 2;
209*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
210*a4bb7517SMichael Lange     snum = fscanf(fd, "%d %d\n", &(ele->nodes[0]), &(ele->nodes[1]));CHKERRQ(snum != ele->numNodes);
21130412aabSMichael Lange     break;
21230412aabSMichael Lange   case 2: /* 3-node triangle */
213*a4bb7517SMichael Lange     ele->dim = 2;
214*a4bb7517SMichael Lange     ele->numNodes = 3;
215*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
216*a4bb7517SMichael Lange     snum = fscanf(fd, "%d %d %d\n", &(ele->nodes[0]), &(ele->nodes[1]), &(ele->nodes[2]));CHKERRQ(snum != ele->numNodes);
21730412aabSMichael Lange     break;
21830412aabSMichael Lange   case 3: /* 4-node quadrangle */
219*a4bb7517SMichael Lange     ele->dim = 2;
220*a4bb7517SMichael Lange     ele->numNodes = 4;
221*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
222*a4bb7517SMichael Lange     snum = fscanf(fd, "%d %d %d %d\n", &(ele->nodes[0]), &(ele->nodes[1]), &(ele->nodes[2]), &(ele->nodes[3]));CHKERRQ(snum != ele->numNodes);
22330412aabSMichael Lange     break;
22430412aabSMichael Lange   case 4: /* 4-node tetrahedron */
225*a4bb7517SMichael Lange     ele->dim  = 3;
226*a4bb7517SMichael Lange     ele->numNodes = 4;
227*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
228*a4bb7517SMichael Lange     snum = fscanf(fd, "%d %d %d %d\n", &(ele->nodes[0]), &(ele->nodes[1]), &(ele->nodes[2]), &(ele->nodes[3]));CHKERRQ(snum != ele->numNodes);
22930412aabSMichael Lange     break;
23030412aabSMichael Lange   case 5: /* 8-node hexahedron */
231*a4bb7517SMichael Lange     ele->dim = 3;
232*a4bb7517SMichael Lange     ele->numNodes = 8;
233*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
234*a4bb7517SMichael Lange     snum = fscanf(fd, "%d %d %d %d %d %d %d %d\n", &(ele->nodes[0]), &(ele->nodes[1]), &(ele->nodes[2]), &(ele->nodes[3]), &(ele->nodes[4]), &(ele->nodes[5]), &(ele->nodes[6]), &(ele->nodes[7]));CHKERRQ(snum != ele->numNodes);
23530412aabSMichael Lange     break;
2366b7f3382SJustin Chang   case 15: /* 1-node vertex */
237*a4bb7517SMichael Lange     ele->dim = 0;
238*a4bb7517SMichael Lange     ele->numNodes = 1;
239*a4bb7517SMichael Lange     ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
240*a4bb7517SMichael Lange     snum = fscanf(fd, "%d\n", &(ele->nodes[0]));CHKERRQ(snum != ele->numNodes);
2416b7f3382SJustin Chang     break;
24230412aabSMichael Lange   default:
24330412aabSMichael Lange     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
24430412aabSMichael Lange   }
24530412aabSMichael Lange   PetscFunctionReturn(0);
246db397164SMichael Lange }
247