xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision f6ac7a6a009568d3077e9d2fef02a9d641ae9d48)
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;
29*f6ac7a6aSMichael Lange   float           version;
307d282ae0SMichael Lange   PetscErrorCode  ierr;
317d282ae0SMichael Lange 
327d282ae0SMichael Lange   PetscFunctionBegin;
33abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
347d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
357d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
367d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
377d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
387d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
39abc86ac4SMichael Lange   if (!rank) {
407d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
417d282ae0SMichael Lange     ierr = PetscViewerRead(vheader, line, 1, PETSC_STRING);CHKERRQ(ierr);
427d282ae0SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
437d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
447d282ae0SMichael Lange     ierr = PetscViewerRead(vheader, line, 2, PETSC_STRING);CHKERRQ(ierr);
45*f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d", &version, &fileType);
46*f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
47*f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
48abc86ac4SMichael Lange   }
49abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
507d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
517d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
527d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
537d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
547d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
557d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
587d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
607d282ae0SMichael Lange   PetscFunctionReturn(0);
617d282ae0SMichael Lange }
627d282ae0SMichael Lange 
637d282ae0SMichael Lange 
647d282ae0SMichael Lange #undef __FUNCT__
65331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
66331830f3SMatthew G. Knepley /*@
677d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
68331830f3SMatthew G. Knepley 
69331830f3SMatthew G. Knepley   Collective on comm
70331830f3SMatthew G. Knepley 
71331830f3SMatthew G. Knepley   Input Parameters:
72331830f3SMatthew G. Knepley + comm  - The MPI communicator
73331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
74331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
75331830f3SMatthew G. Knepley 
76331830f3SMatthew G. Knepley   Output Parameter:
77331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
78331830f3SMatthew G. Knepley 
79331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
80331830f3SMatthew G. Knepley 
81331830f3SMatthew G. Knepley   Level: beginner
82331830f3SMatthew G. Knepley 
83331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
84331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
85331830f3SMatthew G. Knepley @*/
86331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
87331830f3SMatthew G. Knepley {
8804d1ad83SMichael Lange   PetscViewerType vtype;
89a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
90331830f3SMatthew G. Knepley   PetscSection   coordSection;
91331830f3SMatthew G. Knepley   Vec            coordinates;
92331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
93a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
94a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
95331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
96331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
9704d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
98331830f3SMatthew G. Knepley   PetscErrorCode ierr;
99331830f3SMatthew G. Knepley 
100331830f3SMatthew G. Knepley   PetscFunctionBegin;
101331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
102331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
103331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
104331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
10504d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
10604d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
107abc86ac4SMichael Lange   if (!rank || binary) {
108331830f3SMatthew G. Knepley     PetscBool match;
109331830f3SMatthew G. Knepley     int       fileType, dataSize;
110*f6ac7a6aSMichael Lange     float     version;
111331830f3SMatthew G. Knepley 
112331830f3SMatthew G. Knepley     /* Read header */
11304d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
11404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
115331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
11604d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr);
117*f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
118*f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
119*f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
120331830f3SMatthew 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);
12104d1ad83SMichael Lange     if (binary) {
122b9eae255SMichael Lange       int checkInt;
123b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_ENUM);CHKERRQ(ierr);
12404d1ad83SMichael Lange       if (checkInt != 1) {
125b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
12604d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
12704d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
12804d1ad83SMichael Lange       }
1290877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
13004d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
13104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", 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");
133331830f3SMatthew G. Knepley     /* Read vertices */
13404d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
13504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
136331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
13704d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);
1380877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
1390877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
140854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
141331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
142331830f3SMatthew G. Knepley       int i;
143b9eae255SMichael Lange       ierr = PetscViewerRead(viewer, &i, 1, PETSC_ENUM);CHKERRQ(ierr);
144b9eae255SMichael Lange       if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr);
14504d1ad83SMichael Lange       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr);
14604d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
147b9eae255SMichael Lange       if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
148331830f3SMatthew G. Knepley     }
14904d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
15004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", 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");
152331830f3SMatthew G. Knepley     /* Read cells */
15304d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
15404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
155331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
15604d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
1570877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
1580877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
159331830f3SMatthew G. Knepley   }
160db397164SMichael Lange 
161abc86ac4SMichael Lange   if (!rank || binary) {
162a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
163a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
164a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
165a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
166a4bb7517SMichael Lange     ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr);
167a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
16804d1ad83SMichael Lange       ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr);
169a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
170a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
171db397164SMichael Lange     }
17204d1ad83SMichael Lange     ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);;
17304d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
174331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
175331830f3SMatthew G. Knepley   }
176abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
177abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
178a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
179331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
180331830f3SMatthew G. Knepley   if (!rank) {
181a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
182a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
183a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
184a4bb7517SMichael Lange         cell++;
185331830f3SMatthew G. Knepley       }
186331830f3SMatthew G. Knepley     }
187331830f3SMatthew G. Knepley   }
188331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
189a4bb7517SMichael Lange   /* Add cell-vertex connections */
190331830f3SMatthew G. Knepley   if (!rank) {
191331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
192a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
193a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
194a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
195a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
196331830f3SMatthew G. Knepley         }
197a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
198a4bb7517SMichael Lange         cell++;
199331830f3SMatthew G. Knepley       }
200a4bb7517SMichael Lange     }
201331830f3SMatthew G. Knepley   }
2026228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
203c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
204331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
205331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
206331830f3SMatthew G. Knepley   if (interpolate) {
207e56d480eSMatthew G. Knepley     DM idm = NULL;
208331830f3SMatthew G. Knepley 
209331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
210331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
211331830f3SMatthew G. Knepley     *dm  = idm;
212331830f3SMatthew G. Knepley   }
21316ad7e67SMichael Lange 
21416ad7e67SMichael Lange   if (!rank) {
21516ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
21616ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
217d1a54cefSMatthew G. Knepley 
21816ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
21916ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
220a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
22116ad7e67SMichael Lange         PetscInt joinSize;
22216ad7e67SMichael Lange         const PetscInt *join;
223a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
224a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
22516ad7e67SMichael Lange         }
226a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
227a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
228a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
229a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
23016ad7e67SMichael Lange       }
231a4bb7517SMichael Lange     }
23216ad7e67SMichael Lange   }
23316ad7e67SMichael Lange 
234331830f3SMatthew G. Knepley   /* Read coordinates */
235979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
236331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
237331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
23875b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
23975b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
240331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
241331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
242331830f3SMatthew G. Knepley   }
243331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
245331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
246331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
247331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
249331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
250331830f3SMatthew G. Knepley   if (!rank) {
251331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
252331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
253331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
254331830f3SMatthew G. Knepley       }
255331830f3SMatthew G. Knepley     }
256331830f3SMatthew G. Knepley   }
257331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
258331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
259331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
260331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
261a4bb7517SMichael Lange   /* Clean up intermediate storage */
262abc86ac4SMichael Lange   if (!rank || binary) {
263a4bb7517SMichael Lange     for (c = 0; c < numCells; ++c) {
264a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].nodes);
265a4bb7517SMichael Lange       ierr = PetscFree(gmsh_elem[c].tags);
266a4bb7517SMichael Lange     }
267a4bb7517SMichael Lange     ierr = PetscFree(gmsh_elem);
268a4bb7517SMichael Lange   }
269331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
270331830f3SMatthew G. Knepley }
271db397164SMichael Lange 
272db397164SMichael Lange #undef __FUNCT__
27330412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
27404d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele)
275db397164SMichael Lange {
27604d1ad83SMichael Lange   int            cellType, numElem;
277a4bb7517SMichael Lange   PetscErrorCode ierr;
278db397164SMichael Lange 
27930412aabSMichael Lange   PetscFunctionBegin;
28004d1ad83SMichael Lange   if (binary) {
281b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
282b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_ENUM, 1);CHKERRQ(ierr);
283b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_ENUM);CHKERRQ(ierr);
284b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_ENUM, 1);CHKERRQ(ierr);
285b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
286b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_ENUM, 1);CHKERRQ(ierr);
287b9eae255SMichael Lange     ierr = PetscViewerRead(viewer,  &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
288b9eae255SMichael Lange     if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_ENUM, 1);CHKERRQ(ierr);
28904d1ad83SMichael Lange   } else {
290b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_ENUM);CHKERRQ(ierr);
291b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_ENUM);CHKERRQ(ierr);
292b9eae255SMichael Lange     ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_ENUM);CHKERRQ(ierr);
29304d1ad83SMichael Lange   }
294a4bb7517SMichael Lange   ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr);
295b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_ENUM);CHKERRQ(ierr);
296b9eae255SMichael Lange   if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_ENUM, ele->numTags);CHKERRQ(ierr);
29730412aabSMichael Lange   switch (cellType) {
29830412aabSMichael Lange   case 1: /* 2-node line */
299a4bb7517SMichael Lange     ele->dim = 1;
300a4bb7517SMichael Lange     ele->numNodes = 2;
30130412aabSMichael Lange     break;
30230412aabSMichael Lange   case 2: /* 3-node triangle */
303a4bb7517SMichael Lange     ele->dim = 2;
304a4bb7517SMichael Lange     ele->numNodes = 3;
30530412aabSMichael Lange     break;
30630412aabSMichael Lange   case 3: /* 4-node quadrangle */
307a4bb7517SMichael Lange     ele->dim = 2;
308a4bb7517SMichael Lange     ele->numNodes = 4;
30930412aabSMichael Lange     break;
31030412aabSMichael Lange   case 4: /* 4-node tetrahedron */
311a4bb7517SMichael Lange     ele->dim  = 3;
312a4bb7517SMichael Lange     ele->numNodes = 4;
31330412aabSMichael Lange     break;
31430412aabSMichael Lange   case 5: /* 8-node hexahedron */
315a4bb7517SMichael Lange     ele->dim = 3;
316a4bb7517SMichael Lange     ele->numNodes = 8;
31730412aabSMichael Lange     break;
3186b7f3382SJustin Chang   case 15: /* 1-node vertex */
319a4bb7517SMichael Lange     ele->dim = 0;
320a4bb7517SMichael Lange     ele->numNodes = 1;
3216b7f3382SJustin Chang     break;
32230412aabSMichael Lange   default:
32330412aabSMichael Lange     SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
32430412aabSMichael Lange   }
32504d1ad83SMichael Lange   ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr);
326b9eae255SMichael Lange   ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_ENUM);CHKERRQ(ierr);
3270877b519SMichael Lange   if (byteSwap) {ierr = PetscByteSwap(ele->nodes, PETSC_ENUM, ele->numNodes);CHKERRQ(ierr);}
32830412aabSMichael Lange   PetscFunctionReturn(0);
329db397164SMichael Lange }
330