xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision b6a3fe3c56097d5d92c7d822d284fb4cf15ccc65)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #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;
28c1c22fd2SMatthew G. Knepley   int             fT;
29c1c22fd2SMatthew G. Knepley   PetscInt        fileType;
30f6ac7a6aSMichael Lange   float           version;
317d282ae0SMichael Lange   PetscErrorCode  ierr;
327d282ae0SMichael Lange 
337d282ae0SMichael Lange   PetscFunctionBegin;
34abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
357d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
367d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
377d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
387d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
397d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
40abc86ac4SMichael Lange   if (!rank) {
417d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
42060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
437d282ae0SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
447d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
45060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
46c1c22fd2SMatthew G. Knepley     snum = sscanf(line, "%f %d", &version, &fT);
47c1c22fd2SMatthew G. Knepley     fileType = (PetscInt) fT;
48f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
49f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
50abc86ac4SMichael Lange   }
51abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
527d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
537d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
547d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
557d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
587d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
607d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
617d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
627d282ae0SMichael Lange   PetscFunctionReturn(0);
637d282ae0SMichael Lange }
647d282ae0SMichael Lange 
657d282ae0SMichael Lange 
667d282ae0SMichael Lange #undef __FUNCT__
67331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh"
68331830f3SMatthew G. Knepley /*@
697d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
70331830f3SMatthew G. Knepley 
71331830f3SMatthew G. Knepley   Collective on comm
72331830f3SMatthew G. Knepley 
73331830f3SMatthew G. Knepley   Input Parameters:
74331830f3SMatthew G. Knepley + comm  - The MPI communicator
75331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
76331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
77331830f3SMatthew G. Knepley 
78331830f3SMatthew G. Knepley   Output Parameter:
79331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
80331830f3SMatthew G. Knepley 
81331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
823d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
83331830f3SMatthew G. Knepley 
84331830f3SMatthew G. Knepley   Level: beginner
85331830f3SMatthew G. Knepley 
86331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
87331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
88331830f3SMatthew G. Knepley @*/
89331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
90331830f3SMatthew G. Knepley {
9104d1ad83SMichael Lange   PetscViewerType vtype;
92a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
93331830f3SMatthew G. Knepley   PetscSection   coordSection;
94331830f3SMatthew G. Knepley   Vec            coordinates;
95331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
96a4bb7517SMichael Lange   PetscInt       dim = 0, coordSize, c, v, d, cell;
97a4bb7517SMichael Lange   int            numVertices = 0, numCells = 0, trueNumCells = 0, snum;
98331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
99331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
10004d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
101331830f3SMatthew G. Knepley   PetscErrorCode ierr;
102331830f3SMatthew G. Knepley 
103331830f3SMatthew G. Knepley   PetscFunctionBegin;
104331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
105331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
106331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
107331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1083b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
10904d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
11004d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
111abc86ac4SMichael Lange   if (!rank || binary) {
112331830f3SMatthew G. Knepley     PetscBool match;
113331830f3SMatthew G. Knepley     int       fileType, dataSize;
114f6ac7a6aSMichael Lange     float     version;
115331830f3SMatthew G. Knepley 
116331830f3SMatthew G. Knepley     /* Read header */
117060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
11804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
119331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
120060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
121f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
122f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
123f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
124331830f3SMatthew 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);
12504d1ad83SMichael Lange     if (binary) {
126b9eae255SMichael Lange       int checkInt;
127060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
12804d1ad83SMichael Lange       if (checkInt != 1) {
129b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
13004d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
13104d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
13204d1ad83SMichael Lange       }
1330877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
134060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
13504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", 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");
137331830f3SMatthew G. Knepley     /* Read vertices */
138060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
13904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
140331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
141060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
1420877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
1430877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
144854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
145331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
146331830f3SMatthew G. Knepley       int i;
147060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
148b9eae255SMichael Lange       if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr);
149060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
15004d1ad83SMichael Lange       if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr);
151b9eae255SMichael Lange       if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
152331830f3SMatthew G. Knepley     }
153060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
15404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", 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");
156331830f3SMatthew G. Knepley     /* Read cells */
157060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
15804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
159331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
160060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
1610877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
1620877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
163331830f3SMatthew G. Knepley   }
164db397164SMichael Lange 
165abc86ac4SMichael Lange   if (!rank || binary) {
166a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
167a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
168a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
169a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
1703d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
171a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
172a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
173a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
174db397164SMichael Lange     }
175060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
17604d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
177331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
178331830f3SMatthew G. Knepley   }
179abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
180abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
181a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
182331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
183331830f3SMatthew G. Knepley   if (!rank) {
184a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
185a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
186a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
187a4bb7517SMichael Lange         cell++;
188331830f3SMatthew G. Knepley       }
189331830f3SMatthew G. Knepley     }
190331830f3SMatthew G. Knepley   }
191331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
192a4bb7517SMichael Lange   /* Add cell-vertex connections */
193331830f3SMatthew G. Knepley   if (!rank) {
194331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
195a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
196a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
197a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
198a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
199331830f3SMatthew G. Knepley         }
200a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
201a4bb7517SMichael Lange         cell++;
202331830f3SMatthew G. Knepley       }
203a4bb7517SMichael Lange     }
204331830f3SMatthew G. Knepley   }
2056228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
206c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
207331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
208331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
209331830f3SMatthew G. Knepley   if (interpolate) {
210e56d480eSMatthew G. Knepley     DM idm = NULL;
211331830f3SMatthew G. Knepley 
212331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
213331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
214331830f3SMatthew G. Knepley     *dm  = idm;
215331830f3SMatthew G. Knepley   }
21616ad7e67SMichael Lange 
21716ad7e67SMichael Lange   if (!rank) {
21816ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
21916ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
220d1a54cefSMatthew G. Knepley 
22116ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
22216ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
223a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
22416ad7e67SMichael Lange         PetscInt joinSize;
22516ad7e67SMichael Lange         const PetscInt *join;
226a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
227a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
22816ad7e67SMichael Lange         }
229a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
230a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
231a4bb7517SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
232a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
23316ad7e67SMichael Lange       }
234a4bb7517SMichael Lange     }
23516ad7e67SMichael Lange   }
23616ad7e67SMichael Lange 
237331830f3SMatthew G. Knepley   /* Read coordinates */
238979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
239331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
240331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
24175b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
24275b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
243331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
244331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
245331830f3SMatthew G. Knepley   }
246331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
247331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
249331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
250331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
251331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
252331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
253331830f3SMatthew G. Knepley   if (!rank) {
254331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
255331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
256331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
257331830f3SMatthew G. Knepley       }
258331830f3SMatthew G. Knepley     }
259331830f3SMatthew G. Knepley   }
260331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
261331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
262331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
263331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
264a4bb7517SMichael Lange   /* Clean up intermediate storage */
26529403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
2663b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
267331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
268331830f3SMatthew G. Knepley }
269db397164SMichael Lange 
270db397164SMichael Lange #undef __FUNCT__
27130412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
2723d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
273db397164SMichael Lange {
2743d51f2daSMichael Lange   PetscInt       c;
2753d51f2daSMichael Lange   GmshElement   *elements;
276*b6a3fe3cSMichael Lange   int            i, cellType, dim, numNodes, numElem, numTags;
277*b6a3fe3cSMichael Lange   int            ibuf[4];
278a4bb7517SMichael Lange   PetscErrorCode ierr;
279db397164SMichael Lange 
28030412aabSMichael Lange   PetscFunctionBegin;
2813d51f2daSMichael Lange   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
2823d51f2daSMichael Lange   for (c = 0; c < numCells;) {
283*b6a3fe3cSMichael Lange     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
284*b6a3fe3cSMichael Lange     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
28504d1ad83SMichael Lange     if (binary) {
286*b6a3fe3cSMichael Lange       cellType = ibuf[0];
287*b6a3fe3cSMichael Lange       numElem = ibuf[1];
288*b6a3fe3cSMichael Lange       numTags = ibuf[2];
28904d1ad83SMichael Lange     } else {
290*b6a3fe3cSMichael Lange       elements[c].id = ibuf[0];
291*b6a3fe3cSMichael Lange       cellType = ibuf[1];
292*b6a3fe3cSMichael Lange       numTags = ibuf[2];
2933d51f2daSMichael Lange       numElem = 1;
29404d1ad83SMichael Lange     }
295*b6a3fe3cSMichael Lange     switch (cellType) {
296*b6a3fe3cSMichael Lange     case 1: /* 2-node line */
297*b6a3fe3cSMichael Lange       dim = 1;
298*b6a3fe3cSMichael Lange       numNodes = 2;
299*b6a3fe3cSMichael Lange       break;
300*b6a3fe3cSMichael Lange     case 2: /* 3-node triangle */
301*b6a3fe3cSMichael Lange       dim = 2;
302*b6a3fe3cSMichael Lange       numNodes = 3;
303*b6a3fe3cSMichael Lange       break;
304*b6a3fe3cSMichael Lange     case 3: /* 4-node quadrangle */
305*b6a3fe3cSMichael Lange       dim = 2;
306*b6a3fe3cSMichael Lange       numNodes = 4;
307*b6a3fe3cSMichael Lange       break;
308*b6a3fe3cSMichael Lange     case 4: /* 4-node tetrahedron */
309*b6a3fe3cSMichael Lange       dim  = 3;
310*b6a3fe3cSMichael Lange       numNodes = 4;
311*b6a3fe3cSMichael Lange       break;
312*b6a3fe3cSMichael Lange     case 5: /* 8-node hexahedron */
313*b6a3fe3cSMichael Lange       dim = 3;
314*b6a3fe3cSMichael Lange       numNodes = 8;
315*b6a3fe3cSMichael Lange       break;
316*b6a3fe3cSMichael Lange     case 15: /* 1-node vertex */
317*b6a3fe3cSMichael Lange       dim = 0;
318*b6a3fe3cSMichael Lange       numNodes = 1;
319*b6a3fe3cSMichael Lange       break;
320*b6a3fe3cSMichael Lange     default:
321*b6a3fe3cSMichael Lange       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
322*b6a3fe3cSMichael Lange     }
3233d51f2daSMichael Lange     for (i = 0; i < numElem; ++i, ++c) {
3243d51f2daSMichael Lange       /* Loop over inner binary element block */
3253d51f2daSMichael Lange       if (binary) {
3263d51f2daSMichael Lange         ierr = PetscViewerRead(viewer,  &(elements[c].id), 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
3273d51f2daSMichael Lange         if (byteSwap) ierr = PetscByteSwap( &(elements[c].id), PETSC_ENUM, 1);CHKERRQ(ierr);
3283d51f2daSMichael Lange       }
329*b6a3fe3cSMichael Lange       elements[c].dim = dim;
330*b6a3fe3cSMichael Lange       elements[c].numNodes = numNodes;
331*b6a3fe3cSMichael Lange       elements[c].numTags = numTags;
332*b6a3fe3cSMichael Lange 
3333d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
3343d51f2daSMichael Lange       if (byteSwap) ierr = PetscByteSwap(elements[c].tags, PETSC_ENUM, elements[c].numTags);CHKERRQ(ierr);
3353d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
3363d51f2daSMichael Lange       if (byteSwap) {ierr = PetscByteSwap(elements[c].nodes, PETSC_ENUM, elements[c].numNodes);CHKERRQ(ierr);}
3373d51f2daSMichael Lange     }
3383d51f2daSMichael Lange   }
3393d51f2daSMichael Lange   *gmsh_elems = elements;
34030412aabSMichael Lange   PetscFunctionReturn(0);
341db397164SMichael Lange }
342