xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 0877b5199acc9ccee243ec651c7ed0bdd5e0b7ce)
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;
28b9eae255SMichael 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);
44*0877b519SMichael Lange     snum = sscanf(line, "2.2 %d", &fileType);
45*0877b519SMichael Lange     if (snum != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
46abc86ac4SMichael Lange   }
47abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
487d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
497d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
507d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
517d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
527d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
537d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
547d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
557d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
587d282ae0SMichael Lange   PetscFunctionReturn(0);
597d282ae0SMichael Lange }
607d282ae0SMichael Lange 
617d282ae0SMichael Lange 
627d282ae0SMichael Lange #undef __FUNCT__
63331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
64331830f3SMatthew G. Knepley /*@
657d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
66331830f3SMatthew G. Knepley 
67331830f3SMatthew G. Knepley   Collective on comm
68331830f3SMatthew G. Knepley 
69331830f3SMatthew G. Knepley   Input Parameters:
70331830f3SMatthew G. Knepley + comm  - The MPI communicator
71331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
72331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
73331830f3SMatthew G. Knepley 
74331830f3SMatthew G. Knepley   Output Parameter:
75331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
76331830f3SMatthew G. Knepley 
77331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
78331830f3SMatthew G. Knepley 
79331830f3SMatthew G. Knepley   Level: beginner
80331830f3SMatthew G. Knepley 
81331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
82331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
83331830f3SMatthew G. Knepley @*/
84331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
85331830f3SMatthew G. Knepley {
8604d1ad83SMichael Lange   PetscViewerType vtype;
87a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
88331830f3SMatthew G. Knepley   PetscSection   coordSection;
89331830f3SMatthew G. Knepley   Vec            coordinates;
90331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
91a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
92a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
93331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
94331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
9504d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
96331830f3SMatthew G. Knepley   PetscErrorCode ierr;
97331830f3SMatthew G. Knepley 
98331830f3SMatthew G. Knepley   PetscFunctionBegin;
99331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
100331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
101331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
102331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10304d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
10404d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
105abc86ac4SMichael Lange   if (!rank || binary) {
106331830f3SMatthew G. Knepley     PetscBool match;
107331830f3SMatthew G. Knepley     int       fileType, dataSize;
108331830f3SMatthew G. Knepley 
109331830f3SMatthew G. Knepley     /* Read header */
11004d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
11104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
112331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
11304d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr);
114*0877b519SMichael Lange     snum = sscanf(line, "2.2 %d %d", &fileType, &dataSize);
115*0877b519SMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
116331830f3SMatthew 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);
11704d1ad83SMichael Lange     if (binary) {
118b9eae255SMichael Lange       int checkInt;
119b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_ENUM);CHKERRQ(ierr);
12004d1ad83SMichael Lange       if (checkInt != 1) {
121b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
12204d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
12304d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
12404d1ad83SMichael Lange       }
125*0877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
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);
134*0877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
135*0877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
136331830f3SMatthew G. Knepley     ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr);
137331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
138331830f3SMatthew G. Knepley       int i;
139b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &i, 1, PETSC_ENUM);CHKERRQ(ierr);
140b9eae255SMichael Lange       if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr);
14104d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr);
14204d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
143b9eae255SMichael Lange       if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
144331830f3SMatthew G. Knepley     }
14504d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
14604d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
147331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
148331830f3SMatthew G. Knepley     /* Read cells */
14904d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
15004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
151331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
15204d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
153*0877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
154*0877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
155331830f3SMatthew G. Knepley   }
156db397164SMichael Lange 
157abc86ac4SMichael Lange   if (!rank || binary) {
158a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
159a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
160a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
161a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
162a4bb7517SMichael Lange     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
163a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
16404d1ad83SMichael Lange       ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr);
165a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
166a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
167331830f3SMatthew G. Knepley     }
16804d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
16904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
170331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
171331830f3SMatthew G. Knepley   }
172abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
173abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
174a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
175a4bb7517SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
176a4bb7517SMichael Lange   if (!rank) {
177a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
178a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
179a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
180a4bb7517SMichael Lange         cell++;
181a4bb7517SMichael Lange       }
182a4bb7517SMichael Lange     }
183a4bb7517SMichael Lange   }
184a4bb7517SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
185a4bb7517SMichael Lange   /* Add cell-vertex connections */
186a4bb7517SMichael Lange   if (!rank) {
187a4bb7517SMichael Lange     PetscInt pcone[8], corner;
188a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
189a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
190a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
191a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
192a4bb7517SMichael Lange         }
193a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
194a4bb7517SMichael Lange         cell++;
195a4bb7517SMichael Lange       }
196a4bb7517SMichael Lange     }
197a4bb7517SMichael Lange   }
1986228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
199a4bb7517SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
200331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
201331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
202331830f3SMatthew G. Knepley   if (interpolate) {
203e56d480eSMatthew G. Knepley     DM idm = NULL;
204331830f3SMatthew G. Knepley 
205331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
206331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
207331830f3SMatthew G. Knepley     *dm  = idm;
208331830f3SMatthew G. Knepley   }
20916ad7e67SMichael Lange 
21016ad7e67SMichael Lange   if (!rank) {
21116ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
21216ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
213d1a54cefSMatthew G. Knepley 
21416ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
21516ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
216a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
21716ad7e67SMichael Lange         PetscInt joinSize;
21816ad7e67SMichael Lange         const PetscInt *join;
219a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
220a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
22116ad7e67SMichael Lange         }
222a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
223a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
224a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
225a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
22616ad7e67SMichael Lange       }
227a4bb7517SMichael Lange     }
22816ad7e67SMichael Lange   }
22916ad7e67SMichael Lange 
230331830f3SMatthew G. Knepley   /* Read coordinates */
231979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
232331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
233331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
23475b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
23575b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
236331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
237331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
238331830f3SMatthew G. Knepley   }
239331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
240331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
241331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
242331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
243331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
245331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
246331830f3SMatthew G. Knepley   if (!rank) {
247331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
248331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
249331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
250331830f3SMatthew G. Knepley       }
251331830f3SMatthew G. Knepley     }
252331830f3SMatthew G. Knepley   }
253331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
254331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
255331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
256331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
257a4bb7517SMichael Lange   /* Clean up intermediate storage */
258abc86ac4SMichael Lange   if (!rank || binary) {
259a4bb7517SMichael Lange     for (c = 0; c < numCells; ++c) {
260a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].nodes);
261a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].tags);
262a4bb7517SMichael Lange     }
263a4bb7517SMichael Lange     ierr = PetscFree(gmsh_elem);
264a4bb7517SMichael Lange   }
265331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
266331830f3SMatthew G. Knepley }
267db397164SMichael Lange 
268db397164SMichael Lange #undef __FUNCT__
26930412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
27004d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele)
271db397164SMichael Lange {
27204d1ad83SMichael Lange   int            cellType, numElem;
273a4bb7517SMichael Lange   PetscErrorCode ierr;
274db397164SMichael Lange 
27530412aabSMichael Lange   PetscFunctionBegin;
27604d1ad83SMichael Lange   if (binary) {
277b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
278b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_ENUM, 1);CHKERRQ(ierr);
279b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_ENUM);CHKERRQ(ierr);
280b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_ENUM, 1);CHKERRQ(ierr);
281b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
282b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_ENUM, 1);CHKERRQ(ierr);
283b9eae255SMichael Lange     ierr = PetscViewerRead(viewer,  &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
284b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_ENUM, 1);CHKERRQ(ierr);
28504d1ad83SMichael Lange   } else {
286b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
287b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
288b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
28904d1ad83SMichael Lange   }
290a4bb7517SMichael Lange   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
291b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_ENUM);CHKERRQ(ierr);
292b9eae255SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_ENUM, ele->numTags);CHKERRQ(ierr);
29330412aabSMichael Lange   switch (cellType) {
29430412aabSMichael Lange   case 1: /* 2-node line */
295a4bb7517SMichael Lange     ele->dim = 1;
296a4bb7517SMichael Lange     ele->numNodes = 2;
29730412aabSMichael Lange     break;
29830412aabSMichael Lange   case 2: /* 3-node triangle */
299a4bb7517SMichael Lange     ele->dim = 2;
300a4bb7517SMichael Lange     ele->numNodes = 3;
30130412aabSMichael Lange     break;
30230412aabSMichael Lange   case 3: /* 4-node quadrangle */
303a4bb7517SMichael Lange     ele->dim = 2;
304a4bb7517SMichael Lange     ele->numNodes = 4;
30530412aabSMichael Lange     break;
30630412aabSMichael Lange   case 4: /* 4-node tetrahedron */
307a4bb7517SMichael Lange     ele->dim  = 3;
308a4bb7517SMichael Lange     ele->numNodes = 4;
30930412aabSMichael Lange     break;
31030412aabSMichael Lange   case 5: /* 8-node hexahedron */
311a4bb7517SMichael Lange     ele->dim = 3;
312a4bb7517SMichael Lange     ele->numNodes = 8;
31330412aabSMichael Lange     break;
3146b7f3382SJustin Chang   case 15: /* 1-node vertex */
315a4bb7517SMichael Lange     ele->dim = 0;
316a4bb7517SMichael Lange     ele->numNodes = 1;
3176b7f3382SJustin Chang     break;
31830412aabSMichael Lange   default:
31930412aabSMichael Lange     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
32030412aabSMichael Lange   }
32104d1ad83SMichael Lange   ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
322b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_ENUM);CHKERRQ(ierr);
323*0877b519SMichael Lange   if (byteSwap) {ierr = PetscByteSwap(ele->nodes, PETSC_ENUM, ele->numNodes);CHKERRQ(ierr);}
32430412aabSMichael Lange   PetscFunctionReturn(0);
325db397164SMichael Lange }
326