xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 7d282ae044b3c3a4da9af0ebb02c192250bc1e8b)
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__
5*7d282ae0SMichael Lange #define __FUNCT__ "DMPlexCreateGmshFromFile"
6*7d282ae0SMichael Lange /*@C
7*7d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
8*7d282ae0SMichael Lange 
9*7d282ae0SMichael Lange + comm        - The MPI communicator
10*7d282ae0SMichael Lange . filename    - Name of the Gmsh file
11*7d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
12*7d282ae0SMichael Lange 
13*7d282ae0SMichael Lange   Output Parameter:
14*7d282ae0SMichael Lange . dm  - The DM object representing the mesh
15*7d282ae0SMichael Lange 
16*7d282ae0SMichael Lange   Level: beginner
17*7d282ae0SMichael Lange 
18*7d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
19*7d282ae0SMichael Lange @*/
20*7d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
21*7d282ae0SMichael Lange {
22*7d282ae0SMichael Lange   PetscViewer     viewer, vheader;
23*7d282ae0SMichael Lange   PetscViewerType vtype;
24*7d282ae0SMichael Lange   char            line[PETSC_MAX_PATH_LEN];
25*7d282ae0SMichael Lange   int             snum;
26*7d282ae0SMichael Lange   PetscBool       match;
27*7d282ae0SMichael Lange   PetscInt        fileType;
28*7d282ae0SMichael Lange   PetscErrorCode  ierr;
29*7d282ae0SMichael Lange 
30*7d282ae0SMichael Lange   PetscFunctionBegin;
31*7d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
32*7d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
33*7d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
34*7d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
35*7d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
36*7d282ae0SMichael Lange   /* Read only the first two lines of the Gmsh file */
37*7d282ae0SMichael Lange   ierr = PetscViewerRead(vheader, line, 1, PETSC_STRING);CHKERRQ(ierr);
38*7d282ae0SMichael Lange   ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
39*7d282ae0SMichael Lange   if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
40*7d282ae0SMichael Lange   ierr = PetscViewerRead(vheader, line, 2, PETSC_STRING);CHKERRQ(ierr);
41*7d282ae0SMichael Lange   snum = sscanf(line, "2.2 %d", &fileType);CHKERRQ(snum != 1);
42*7d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
43*7d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
44*7d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
45*7d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
46*7d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
47*7d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
48*7d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
49*7d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
50*7d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
51*7d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
52*7d282ae0SMichael Lange   PetscFunctionReturn(0);
53*7d282ae0SMichael Lange }
54*7d282ae0SMichael Lange 
55*7d282ae0SMichael Lange 
56*7d282ae0SMichael Lange #undef __FUNCT__
57331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
58331830f3SMatthew G. Knepley /*@
59*7d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
60331830f3SMatthew G. Knepley 
61331830f3SMatthew G. Knepley   Collective on comm
62331830f3SMatthew G. Knepley 
63331830f3SMatthew G. Knepley   Input Parameters:
64331830f3SMatthew G. Knepley + comm  - The MPI communicator
65331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
66331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
67331830f3SMatthew G. Knepley 
68331830f3SMatthew G. Knepley   Output Parameter:
69331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
70331830f3SMatthew G. Knepley 
71331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
72331830f3SMatthew G. Knepley 
73331830f3SMatthew G. Knepley   Level: beginner
74331830f3SMatthew G. Knepley 
75331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
76331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
77331830f3SMatthew G. Knepley @*/
78331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
79331830f3SMatthew G. Knepley {
8004d1ad83SMichael Lange   PetscViewerType vtype;
81a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
82331830f3SMatthew G. Knepley   PetscSection   coordSection;
83331830f3SMatthew G. Knepley   Vec            coordinates;
84331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
85a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
86a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
87331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
88331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
8904d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
90331830f3SMatthew G. Knepley   PetscErrorCode ierr;
91331830f3SMatthew G. Knepley 
92331830f3SMatthew G. Knepley   PetscFunctionBegin;
93331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
94331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
95331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
96331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
9704d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
9804d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
99331830f3SMatthew G. Knepley   if (!rank) {
100331830f3SMatthew G. Knepley     PetscBool match;
101331830f3SMatthew G. Knepley     int       fileType, dataSize;
102331830f3SMatthew G. Knepley 
103331830f3SMatthew G. Knepley     /* Read header */
10404d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
10504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
106331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
10704d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr);
10804d1ad83SMichael Lange     snum = sscanf(line, "2.2 %d %d", &fileType, &dataSize);CHKERRQ(snum != 2);
109331830f3SMatthew 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);
11004d1ad83SMichael Lange     if (binary) {
11104d1ad83SMichael Lange       PetscInt checkInt;
11204d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_INT);CHKERRQ(ierr);
11304d1ad83SMichael Lange       if (checkInt != 1) {
11404d1ad83SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_INT, 1);CHKERRQ(ierr);
11504d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
11604d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
11704d1ad83SMichael Lange       }
11804d1ad83SMichael Lange     } else {
11904d1ad83SMichael Lange       if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
12004d1ad83SMichael Lange     }
12104d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
12204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
123331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
124331830f3SMatthew G. Knepley     /* Read vertices */
12504d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
12604d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
127331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
12804d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
12904d1ad83SMichael Lange     snum = sscanf(line, "%d", &numVertices);CHKERRQ(snum != 1);
130331830f3SMatthew G. Knepley     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
131331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
132331830f3SMatthew G. Knepley       int    i;
13304d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &i, 1, PETSC_INT);CHKERRQ(ierr);
13404d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&i, PETSC_INT, 1);CHKERRQ(ierr);
13504d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr);
13604d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
137331830f3SMatthew G. Knepley       if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
138331830f3SMatthew G. Knepley     }
13904d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
141331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
142331830f3SMatthew G. Knepley     /* Read cells */
14304d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
145331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
14604d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14704d1ad83SMichael Lange     snum = sscanf(line, "%d", &numCells);CHKERRQ(snum != 1);
148331830f3SMatthew G. Knepley   }
149db397164SMichael Lange 
150331830f3SMatthew G. Knepley   if (!rank) {
151a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
152a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
153a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
154a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
155a4bb7517SMichael Lange     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
156a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
15704d1ad83SMichael Lange       ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr);
158a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
159a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
160331830f3SMatthew G. Knepley     }
16104d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
16204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
163331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
164331830f3SMatthew G. Knepley   }
165a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
166a4bb7517SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
167a4bb7517SMichael Lange   if (!rank) {
168a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
169a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
170a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
171a4bb7517SMichael Lange         cell++;
172a4bb7517SMichael Lange       }
173a4bb7517SMichael Lange     }
174a4bb7517SMichael Lange   }
175a4bb7517SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
176a4bb7517SMichael Lange   /* Add cell-vertex connections */
177a4bb7517SMichael Lange   if (!rank) {
178a4bb7517SMichael Lange     PetscInt pcone[8], corner;
179a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
180a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
181a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
182a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
183a4bb7517SMichael Lange         }
184a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
185a4bb7517SMichael Lange         cell++;
186a4bb7517SMichael Lange       }
187a4bb7517SMichael Lange     }
188a4bb7517SMichael Lange   }
1896228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
190a4bb7517SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
191331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
192331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
193331830f3SMatthew G. Knepley   if (interpolate) {
194e56d480eSMatthew G. Knepley     DM idm = NULL;
195331830f3SMatthew G. Knepley 
196331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
197331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
198331830f3SMatthew G. Knepley     *dm  = idm;
199331830f3SMatthew G. Knepley   }
20016ad7e67SMichael Lange 
20116ad7e67SMichael Lange   if (!rank) {
20216ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
20316ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
204d1a54cefSMatthew G. Knepley 
20516ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
20616ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
207a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
20816ad7e67SMichael Lange         PetscInt joinSize;
20916ad7e67SMichael Lange         const PetscInt *join;
210a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
211a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
21216ad7e67SMichael Lange         }
213a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
214a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
215a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
216a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
21716ad7e67SMichael Lange       }
218a4bb7517SMichael Lange     }
21916ad7e67SMichael Lange   }
22016ad7e67SMichael Lange 
221331830f3SMatthew G. Knepley   /* Read coordinates */
222979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
223331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
224331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
22575b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
22675b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
227331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
228331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
229331830f3SMatthew G. Knepley   }
230331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
231331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
232331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
233331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
234331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
235331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
236331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
237331830f3SMatthew G. Knepley   if (!rank) {
238331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
239331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
240331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
241331830f3SMatthew G. Knepley       }
242331830f3SMatthew G. Knepley     }
243331830f3SMatthew G. Knepley   }
244331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
245331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
246331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
247331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
248a4bb7517SMichael Lange   /* Clean up intermediate storage */
249a4bb7517SMichael Lange   if (!rank) {
250a4bb7517SMichael Lange     for (c = 0; c < numCells; ++c) {
251a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].nodes);
252a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].tags);
253a4bb7517SMichael Lange     }
254a4bb7517SMichael Lange     ierr = PetscFree(gmsh_elem);
255a4bb7517SMichael Lange   }
256331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
257331830f3SMatthew G. Knepley }
258db397164SMichael Lange 
259db397164SMichael Lange #undef __FUNCT__
26030412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
26104d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele)
262db397164SMichael Lange {
26304d1ad83SMichael Lange   int            cellType, numElem;
264a4bb7517SMichael Lange   PetscErrorCode ierr;
265db397164SMichael Lange 
26630412aabSMichael Lange   PetscFunctionBegin;
26704d1ad83SMichael Lange   if (binary) {
26804d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_INT);CHKERRQ(ierr);
26904d1ad83SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_INT, 1);CHKERRQ(ierr);
27004d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_INT);CHKERRQ(ierr);
27104d1ad83SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_INT, 1);CHKERRQ(ierr);
27204d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_INT);CHKERRQ(ierr);
27304d1ad83SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_INT, 1);CHKERRQ(ierr);
27404d1ad83SMichael Lange     ierr = PetscViewerRead(viewer,  &(ele->id), 1, PETSC_INT);CHKERRQ(ierr);
27504d1ad83SMichael Lange     if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_INT, 1);CHKERRQ(ierr);
27604d1ad83SMichael Lange   } else {
27704d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_INT);CHKERRQ(ierr);
27804d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_INT);CHKERRQ(ierr);
27904d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_INT);CHKERRQ(ierr);
28004d1ad83SMichael Lange   }
281a4bb7517SMichael Lange   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
28204d1ad83SMichael Lange   ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_INT);CHKERRQ(ierr);
28304d1ad83SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_INT, ele->numTags);CHKERRQ(ierr);
28430412aabSMichael Lange   switch (cellType) {
28530412aabSMichael Lange   case 1: /* 2-node line */
286a4bb7517SMichael Lange     ele->dim = 1;
287a4bb7517SMichael Lange     ele->numNodes = 2;
28830412aabSMichael Lange     break;
28930412aabSMichael Lange   case 2: /* 3-node triangle */
290a4bb7517SMichael Lange     ele->dim = 2;
291a4bb7517SMichael Lange     ele->numNodes = 3;
29230412aabSMichael Lange     break;
29330412aabSMichael Lange   case 3: /* 4-node quadrangle */
294a4bb7517SMichael Lange     ele->dim = 2;
295a4bb7517SMichael Lange     ele->numNodes = 4;
29630412aabSMichael Lange     break;
29730412aabSMichael Lange   case 4: /* 4-node tetrahedron */
298a4bb7517SMichael Lange     ele->dim  = 3;
299a4bb7517SMichael Lange     ele->numNodes = 4;
30030412aabSMichael Lange     break;
30130412aabSMichael Lange   case 5: /* 8-node hexahedron */
302a4bb7517SMichael Lange     ele->dim = 3;
303a4bb7517SMichael Lange     ele->numNodes = 8;
30430412aabSMichael Lange     break;
3056b7f3382SJustin Chang   case 15: /* 1-node vertex */
306a4bb7517SMichael Lange     ele->dim = 0;
307a4bb7517SMichael Lange     ele->numNodes = 1;
3086b7f3382SJustin Chang     break;
30930412aabSMichael Lange   default:
31030412aabSMichael Lange     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
31130412aabSMichael Lange   }
31204d1ad83SMichael Lange   ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
31304d1ad83SMichael Lange   ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_INT);CHKERRQ(ierr);
31404d1ad83SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->nodes, PETSC_INT, ele->numNodes);CHKERRQ(ierr);
31530412aabSMichael Lange   PetscFunctionReturn(0);
316db397164SMichael Lange }
317