xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision b9eae255e52efddd1c632647130783c4af38cda0)
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__
57d282ae0SMichael Lange #define __FUNCT__ "DMPlexCreateGmshFromFile"
67d282ae0SMichael Lange /*@C
77d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
87d282ae0SMichael Lange 
97d282ae0SMichael Lange + comm        - The MPI communicator
107d282ae0SMichael Lange . filename    - Name of the Gmsh file
117d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
127d282ae0SMichael Lange 
137d282ae0SMichael Lange   Output Parameter:
147d282ae0SMichael Lange . dm  - The DM object representing the mesh
157d282ae0SMichael Lange 
167d282ae0SMichael Lange   Level: beginner
177d282ae0SMichael Lange 
187d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
197d282ae0SMichael Lange @*/
207d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
217d282ae0SMichael Lange {
227d282ae0SMichael Lange   PetscViewer     viewer, vheader;
23abc86ac4SMichael Lange   PetscMPIInt     rank;
247d282ae0SMichael Lange   PetscViewerType vtype;
257d282ae0SMichael Lange   char            line[PETSC_MAX_PATH_LEN];
267d282ae0SMichael Lange   int             snum;
277d282ae0SMichael Lange   PetscBool       match;
28*b9eae255SMichael Lange   int             fileType;
297d282ae0SMichael Lange   PetscErrorCode  ierr;
307d282ae0SMichael Lange 
317d282ae0SMichael Lange   PetscFunctionBegin;
32abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
337d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
347d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
357d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
367d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
377d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
38abc86ac4SMichael Lange   if (!rank) {
397d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
407d282ae0SMichael Lange     ierr = PetscViewerRead(vheader, line, 1, PETSC_STRING);CHKERRQ(ierr);
417d282ae0SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
427d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
437d282ae0SMichael Lange     ierr = PetscViewerRead(vheader, line, 2, PETSC_STRING);CHKERRQ(ierr);
447d282ae0SMichael Lange     snum = sscanf(line, "2.2 %d", &fileType);CHKERRQ(snum != 1);
45abc86ac4SMichael Lange   }
46abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
477d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
487d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
497d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
507d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
517d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
527d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
537d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
547d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
557d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
577d282ae0SMichael Lange   PetscFunctionReturn(0);
587d282ae0SMichael Lange }
597d282ae0SMichael Lange 
607d282ae0SMichael Lange 
617d282ae0SMichael Lange #undef __FUNCT__
62331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
63331830f3SMatthew G. Knepley /*@
647d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
65331830f3SMatthew G. Knepley 
66331830f3SMatthew G. Knepley   Collective on comm
67331830f3SMatthew G. Knepley 
68331830f3SMatthew G. Knepley   Input Parameters:
69331830f3SMatthew G. Knepley + comm  - The MPI communicator
70331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
71331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
72331830f3SMatthew G. Knepley 
73331830f3SMatthew G. Knepley   Output Parameter:
74331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
75331830f3SMatthew G. Knepley 
76331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
77331830f3SMatthew G. Knepley 
78331830f3SMatthew G. Knepley   Level: beginner
79331830f3SMatthew G. Knepley 
80331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
81331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
82331830f3SMatthew G. Knepley @*/
83331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
84331830f3SMatthew G. Knepley {
8504d1ad83SMichael Lange   PetscViewerType vtype;
86a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
87331830f3SMatthew G. Knepley   PetscSection   coordSection;
88331830f3SMatthew G. Knepley   Vec            coordinates;
89331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
90a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
91a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
92331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
93331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
9404d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
95331830f3SMatthew G. Knepley   PetscErrorCode ierr;
96331830f3SMatthew G. Knepley 
97331830f3SMatthew G. Knepley   PetscFunctionBegin;
98331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
99331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
100331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
101331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10204d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
10304d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
104abc86ac4SMichael Lange   if (!rank || binary) {
105331830f3SMatthew G. Knepley     PetscBool match;
106331830f3SMatthew G. Knepley     int       fileType, dataSize;
107331830f3SMatthew G. Knepley 
108331830f3SMatthew G. Knepley     /* Read header */
10904d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
11004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
111331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
11204d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr);
11304d1ad83SMichael Lange     snum = sscanf(line, "2.2 %d %d", &fileType, &dataSize);CHKERRQ(snum != 2);
114331830f3SMatthew 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);
11504d1ad83SMichael Lange     if (binary) {
116*b9eae255SMichael Lange       int checkInt;
117*b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_ENUM);CHKERRQ(ierr);
11804d1ad83SMichael Lange       if (checkInt != 1) {
119*b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
12004d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
12104d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
12204d1ad83SMichael Lange       }
12304d1ad83SMichael Lange     } else {
12404d1ad83SMichael Lange       if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
12504d1ad83SMichael Lange     }
12604d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
12704d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
128331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
129331830f3SMatthew G. Knepley     /* Read vertices */
13004d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
13104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
132331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
13304d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
13404d1ad83SMichael Lange     snum = sscanf(line, "%d", &numVertices);CHKERRQ(snum != 1);
135331830f3SMatthew G. Knepley     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
136331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
137331830f3SMatthew G. Knepley       int i;
138*b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &i, 1, PETSC_ENUM);CHKERRQ(ierr);
139*b9eae255SMichael Lange       if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr);
14004d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr);
14104d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
142*b9eae255SMichael Lange       if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
143331830f3SMatthew G. Knepley     }
14404d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
146331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
147331830f3SMatthew G. Knepley     /* Read cells */
14804d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
150331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
15104d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
15204d1ad83SMichael Lange     snum = sscanf(line, "%d", &numCells);CHKERRQ(snum != 1);
153331830f3SMatthew G. Knepley   }
154db397164SMichael Lange 
155abc86ac4SMichael Lange   if (!rank || binary) {
156a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
157a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
158a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
159a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
160a4bb7517SMichael Lange     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
161a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
16204d1ad83SMichael Lange       ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr);
163a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
164a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
165331830f3SMatthew G. Knepley     }
16604d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
16704d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
168331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
169331830f3SMatthew G. Knepley   }
170abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
171abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
172a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
173a4bb7517SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
174a4bb7517SMichael Lange   if (!rank) {
175a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
176a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
177a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
178a4bb7517SMichael Lange         cell++;
179a4bb7517SMichael Lange       }
180a4bb7517SMichael Lange     }
181a4bb7517SMichael Lange   }
182a4bb7517SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
183a4bb7517SMichael Lange   /* Add cell-vertex connections */
184a4bb7517SMichael Lange   if (!rank) {
185a4bb7517SMichael Lange     PetscInt pcone[8], corner;
186a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
187a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
188a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
189a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
190a4bb7517SMichael Lange         }
191a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
192a4bb7517SMichael Lange         cell++;
193a4bb7517SMichael Lange       }
194a4bb7517SMichael Lange     }
195a4bb7517SMichael Lange   }
1966228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
197a4bb7517SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
198331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
199331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
200331830f3SMatthew G. Knepley   if (interpolate) {
201e56d480eSMatthew G. Knepley     DM idm = NULL;
202331830f3SMatthew G. Knepley 
203331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
204331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
205331830f3SMatthew G. Knepley     *dm  = idm;
206331830f3SMatthew G. Knepley   }
20716ad7e67SMichael Lange 
20816ad7e67SMichael Lange   if (!rank) {
20916ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
21016ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
211d1a54cefSMatthew G. Knepley 
21216ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
21316ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
214a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
21516ad7e67SMichael Lange         PetscInt joinSize;
21616ad7e67SMichael Lange         const PetscInt *join;
217a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
218a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
21916ad7e67SMichael Lange         }
220a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
221a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
222a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
223a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
22416ad7e67SMichael Lange       }
225a4bb7517SMichael Lange     }
22616ad7e67SMichael Lange   }
22716ad7e67SMichael Lange 
228331830f3SMatthew G. Knepley   /* Read coordinates */
229979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
230331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
231331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
23275b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
23375b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
234331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
235331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
236331830f3SMatthew G. Knepley   }
237331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
238331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
239331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
240331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
241331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
242331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
243331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley   if (!rank) {
245331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
246331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
247331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
248331830f3SMatthew G. Knepley       }
249331830f3SMatthew G. Knepley     }
250331830f3SMatthew G. Knepley   }
251331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
252331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
253331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
254331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
255a4bb7517SMichael Lange   /* Clean up intermediate storage */
256abc86ac4SMichael Lange   if (!rank || binary) {
257a4bb7517SMichael Lange     for (c = 0; c < numCells; ++c) {
258a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].nodes);
259a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].tags);
260a4bb7517SMichael Lange     }
261a4bb7517SMichael Lange     ierr = PetscFree(gmsh_elem);
262a4bb7517SMichael Lange   }
263331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
264331830f3SMatthew G. Knepley }
265db397164SMichael Lange 
266db397164SMichael Lange #undef __FUNCT__
26730412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
26804d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele)
269db397164SMichael Lange {
27004d1ad83SMichael Lange   int            cellType, numElem;
271a4bb7517SMichael Lange   PetscErrorCode ierr;
272db397164SMichael Lange 
27330412aabSMichael Lange   PetscFunctionBegin;
27404d1ad83SMichael Lange   if (binary) {
275*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
276*b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_ENUM, 1);CHKERRQ(ierr);
277*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_ENUM);CHKERRQ(ierr);
278*b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_ENUM, 1);CHKERRQ(ierr);
279*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
280*b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_ENUM, 1);CHKERRQ(ierr);
281*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer,  &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
282*b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_ENUM, 1);CHKERRQ(ierr);
28304d1ad83SMichael Lange   } else {
284*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
285*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
286*b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
28704d1ad83SMichael Lange   }
288a4bb7517SMichael Lange   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
289*b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_ENUM);CHKERRQ(ierr);
290*b9eae255SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_ENUM, ele->numTags);CHKERRQ(ierr);
29130412aabSMichael Lange   switch (cellType) {
29230412aabSMichael Lange   case 1: /* 2-node line */
293a4bb7517SMichael Lange     ele->dim = 1;
294a4bb7517SMichael Lange     ele->numNodes = 2;
29530412aabSMichael Lange     break;
29630412aabSMichael Lange   case 2: /* 3-node triangle */
297a4bb7517SMichael Lange     ele->dim = 2;
298a4bb7517SMichael Lange     ele->numNodes = 3;
29930412aabSMichael Lange     break;
30030412aabSMichael Lange   case 3: /* 4-node quadrangle */
301a4bb7517SMichael Lange     ele->dim = 2;
302a4bb7517SMichael Lange     ele->numNodes = 4;
30330412aabSMichael Lange     break;
30430412aabSMichael Lange   case 4: /* 4-node tetrahedron */
305a4bb7517SMichael Lange     ele->dim  = 3;
306a4bb7517SMichael Lange     ele->numNodes = 4;
30730412aabSMichael Lange     break;
30830412aabSMichael Lange   case 5: /* 8-node hexahedron */
309a4bb7517SMichael Lange     ele->dim = 3;
310a4bb7517SMichael Lange     ele->numNodes = 8;
31130412aabSMichael Lange     break;
3126b7f3382SJustin Chang   case 15: /* 1-node vertex */
313a4bb7517SMichael Lange     ele->dim = 0;
314a4bb7517SMichael Lange     ele->numNodes = 1;
3156b7f3382SJustin Chang     break;
31630412aabSMichael Lange   default:
31730412aabSMichael Lange     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
31830412aabSMichael Lange   }
31904d1ad83SMichael Lange   ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
320*b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_ENUM);CHKERRQ(ierr);
321*b9eae255SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->nodes, PETSC_ENUM, ele->numNodes);CHKERRQ(ierr);
32230412aabSMichael Lange   PetscFunctionReturn(0);
323db397164SMichael Lange }
324