xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision bf6ba3a332b336b825bc993a29c1acb8ada8d525)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3331830f3SMatthew G. Knepley 
47d282ae0SMichael Lange /*@C
57d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
67d282ae0SMichael Lange 
77d282ae0SMichael Lange + comm        - The MPI communicator
87d282ae0SMichael Lange . filename    - Name of the Gmsh file
97d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
107d282ae0SMichael Lange 
117d282ae0SMichael Lange   Output Parameter:
127d282ae0SMichael Lange . dm  - The DM object representing the mesh
137d282ae0SMichael Lange 
147d282ae0SMichael Lange   Level: beginner
157d282ae0SMichael Lange 
167d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
177d282ae0SMichael Lange @*/
187d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
197d282ae0SMichael Lange {
207d282ae0SMichael Lange   PetscViewer     viewer, vheader;
21abc86ac4SMichael Lange   PetscMPIInt     rank;
227d282ae0SMichael Lange   PetscViewerType vtype;
237d282ae0SMichael Lange   char            line[PETSC_MAX_PATH_LEN];
247d282ae0SMichael Lange   int             snum;
257d282ae0SMichael Lange   PetscBool       match;
26c1c22fd2SMatthew G. Knepley   int             fT;
27c1c22fd2SMatthew G. Knepley   PetscInt        fileType;
28f6ac7a6aSMichael Lange   float           version;
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 */
40060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, 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");
43060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
44c1c22fd2SMatthew G. Knepley     snum = sscanf(line, "%f %d", &version, &fT);
45c1c22fd2SMatthew G. Knepley     fileType = (PetscInt) fT;
46f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
47f6ac7a6aSMichael 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 
63df032b82SMatthew G. Knepley static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
64df032b82SMatthew G. Knepley {
65df032b82SMatthew G. Knepley   PetscInt       c, p;
66df032b82SMatthew G. Knepley   GmshElement   *elements;
67*bf6ba3a3SMatthew G. Knepley   int            i, cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
68*bf6ba3a3SMatthew G. Knepley   PetscInt       pibuf[64];
69df032b82SMatthew G. Knepley   int            ibuf[16];
70df032b82SMatthew G. Knepley   PetscErrorCode ierr;
71df032b82SMatthew G. Knepley 
72df032b82SMatthew G. Knepley   PetscFunctionBegin;
73df032b82SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
74df032b82SMatthew G. Knepley   for (c = 0; c < numCells;) {
75df032b82SMatthew G. Knepley     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
76df032b82SMatthew G. Knepley     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
77df032b82SMatthew G. Knepley     if (binary) {
78df032b82SMatthew G. Knepley       cellType = ibuf[0];
79df032b82SMatthew G. Knepley       numElem = ibuf[1];
80df032b82SMatthew G. Knepley       numTags = ibuf[2];
81df032b82SMatthew G. Knepley     } else {
82df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
83df032b82SMatthew G. Knepley       cellType = ibuf[1];
84df032b82SMatthew G. Knepley       numTags = ibuf[2];
85df032b82SMatthew G. Knepley       numElem = 1;
86df032b82SMatthew G. Knepley     }
87*bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
88*bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
89df032b82SMatthew G. Knepley     switch (cellType) {
90df032b82SMatthew G. Knepley     case 1: /* 2-node line */
91df032b82SMatthew G. Knepley       dim = 1;
92df032b82SMatthew G. Knepley       numNodes = 2;
93df032b82SMatthew G. Knepley       break;
94df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
95df032b82SMatthew G. Knepley       dim = 2;
96df032b82SMatthew G. Knepley       numNodes = 3;
97df032b82SMatthew G. Knepley       break;
98df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
99df032b82SMatthew G. Knepley       dim = 2;
100df032b82SMatthew G. Knepley       numNodes = 4;
101df032b82SMatthew G. Knepley       break;
102df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
103df032b82SMatthew G. Knepley       dim  = 3;
104df032b82SMatthew G. Knepley       numNodes = 4;
105df032b82SMatthew G. Knepley       break;
106df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
107df032b82SMatthew G. Knepley       dim = 3;
108df032b82SMatthew G. Knepley       numNodes = 8;
109df032b82SMatthew G. Knepley       break;
110*bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
111*bf6ba3a3SMatthew G. Knepley       dim = 1;
112*bf6ba3a3SMatthew G. Knepley       numNodes = 2;
113*bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
114*bf6ba3a3SMatthew G. Knepley       break;
115*bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
116*bf6ba3a3SMatthew G. Knepley       dim = 2;
117*bf6ba3a3SMatthew G. Knepley       numNodes = 3;
118*bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
119*bf6ba3a3SMatthew G. Knepley       break;
120df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
121df032b82SMatthew G. Knepley       dim = 0;
122df032b82SMatthew G. Knepley       numNodes = 1;
123df032b82SMatthew G. Knepley       break;
124*bf6ba3a3SMatthew G. Knepley     case 6: /* 6-node prism */
125*bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
126*bf6ba3a3SMatthew G. Knepley     case 10: /* 9-node 2nd order quadrangle */
127*bf6ba3a3SMatthew G. Knepley     case 11: /* 10-node 2nd order tetrahedron */
128*bf6ba3a3SMatthew G. Knepley     case 12: /* 27-node 2nd order hexhedron */
129*bf6ba3a3SMatthew G. Knepley     case 13: /* 19-node 2nd order prism */
130*bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
131df032b82SMatthew G. Knepley     default:
132df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
133df032b82SMatthew G. Knepley     }
134df032b82SMatthew G. Knepley     if (binary) {
135*bf6ba3a3SMatthew G. Knepley       const PetscInt nint = numNodes + numTags + 1 + numNodesIgnore;
136df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
137df032b82SMatthew G. Knepley         /* Loop over inner binary element block */
138df032b82SMatthew G. Knepley         elements[c].dim = dim;
139df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
140df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
141df032b82SMatthew G. Knepley 
142df032b82SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
143df032b82SMatthew G. Knepley         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
144df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
145df032b82SMatthew G. Knepley         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
146df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
147df032b82SMatthew G. Knepley       }
148df032b82SMatthew G. Knepley     } else {
149df032b82SMatthew G. Knepley       elements[c].dim = dim;
150df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
151df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
152df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
153df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
154*bf6ba3a3SMatthew G. Knepley       ierr = PetscViewerRead(viewer, pibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
155df032b82SMatthew G. Knepley       c++;
156df032b82SMatthew G. Knepley     }
157df032b82SMatthew G. Knepley   }
158df032b82SMatthew G. Knepley   *gmsh_elems = elements;
159df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
160df032b82SMatthew G. Knepley }
1617d282ae0SMichael Lange 
162331830f3SMatthew G. Knepley /*@
1637d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
164331830f3SMatthew G. Knepley 
165331830f3SMatthew G. Knepley   Collective on comm
166331830f3SMatthew G. Knepley 
167331830f3SMatthew G. Knepley   Input Parameters:
168331830f3SMatthew G. Knepley + comm  - The MPI communicator
169331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
170331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
171331830f3SMatthew G. Knepley 
172331830f3SMatthew G. Knepley   Output Parameter:
173331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
174331830f3SMatthew G. Knepley 
175331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
1763d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
177331830f3SMatthew G. Knepley 
178331830f3SMatthew G. Knepley   Level: beginner
179331830f3SMatthew G. Knepley 
180331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
181331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
182331830f3SMatthew G. Knepley @*/
183331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
184331830f3SMatthew G. Knepley {
18504d1ad83SMichael Lange   PetscViewerType vtype;
1863c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
187331830f3SMatthew G. Knepley   PetscSection   coordSection;
188331830f3SMatthew G. Knepley   Vec            coordinates;
189331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
190dc0ede3bSMatthew G. Knepley   PetscInt       dim = 0, coordSize, c, v, d, r, cell;
191dc0ede3bSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
192331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
193331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
19404d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
195331830f3SMatthew G. Knepley   PetscErrorCode ierr;
196331830f3SMatthew G. Knepley 
197331830f3SMatthew G. Knepley   PetscFunctionBegin;
198331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
199331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
200331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
201331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2023b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
20304d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
20404d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
205abc86ac4SMichael Lange   if (!rank || binary) {
206331830f3SMatthew G. Knepley     PetscBool match;
207331830f3SMatthew G. Knepley     int       fileType, dataSize;
208f6ac7a6aSMichael Lange     float     version;
209331830f3SMatthew G. Knepley 
210331830f3SMatthew G. Knepley     /* Read header */
211060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
21204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
213331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
214060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
215f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
216f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
217f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
218331830f3SMatthew 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);
21904d1ad83SMichael Lange     if (binary) {
220b9eae255SMichael Lange       int checkInt;
221060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
22204d1ad83SMichael Lange       if (checkInt != 1) {
223b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
22404d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
22504d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
22604d1ad83SMichael Lange       }
2270877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
228060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
22904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
230331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
231dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
232060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
233dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
234dc0ede3bSMatthew G. Knepley     if (match) {
235dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
236dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
237dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
238a8667449SSander Arens       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
239dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
240dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
241dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
242dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
243dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
244dc0ede3bSMatthew G. Knepley     }
245dc0ede3bSMatthew G. Knepley     /* Read vertices */
24604d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
247331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
248060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
2490877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
2500877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
251854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
25213aecfe5SMichael Lange     if (binary) {
25313aecfe5SMichael Lange       size_t doubleSize, intSize;
25413aecfe5SMichael Lange       PetscInt elementSize;
25513aecfe5SMichael Lange       char *buffer;
25613aecfe5SMichael Lange       PetscScalar *baseptr;
25713aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
25813aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
25913aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
26013aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
26113aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
26213aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
263331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
26413aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
26513aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
26613aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
26713aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
26813aecfe5SMichael Lange       }
26913aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
27013aecfe5SMichael Lange     } else {
27113aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
272060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
273060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
274b9eae255SMichael Lange         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
275331830f3SMatthew G. Knepley       }
27613aecfe5SMichael Lange     }
277060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
27804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
279331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
280331830f3SMatthew G. Knepley     /* Read cells */
281060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
283331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
284060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2850877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
2860877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
287331830f3SMatthew G. Knepley   }
288db397164SMichael Lange 
289abc86ac4SMichael Lange   if (!rank || binary) {
290a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
291a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
292a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
293a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
2943d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
295a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
296a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
297a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
298db397164SMichael Lange     }
299060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
3001b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
301331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
302331830f3SMatthew G. Knepley   }
303abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
304abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
305a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
306331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
307331830f3SMatthew G. Knepley   if (!rank) {
308a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
309a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
310a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
311a4bb7517SMichael Lange         cell++;
312331830f3SMatthew G. Knepley       }
313331830f3SMatthew G. Knepley     }
314331830f3SMatthew G. Knepley   }
315331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
316a4bb7517SMichael Lange   /* Add cell-vertex connections */
317331830f3SMatthew G. Knepley   if (!rank) {
318331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
319a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
320a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
321a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
322a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
323331830f3SMatthew G. Knepley         }
32497ed6be6Ssarens         if (dim == 3) {
32597ed6be6Ssarens           /* Tetrahedra are inverted */
32697ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
32797ed6be6Ssarens             PetscInt tmp = pcone[0];
32897ed6be6Ssarens             pcone[0] = pcone[1];
32997ed6be6Ssarens             pcone[1] = tmp;
33097ed6be6Ssarens           }
33197ed6be6Ssarens           /* Hexahedra are inverted */
33297ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
33397ed6be6Ssarens             PetscInt tmp = pcone[1];
33497ed6be6Ssarens             pcone[1] = pcone[3];
33597ed6be6Ssarens             pcone[3] = tmp;
33697ed6be6Ssarens           }
33797ed6be6Ssarens         }
338a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
339a4bb7517SMichael Lange         cell++;
340331830f3SMatthew G. Knepley       }
341a4bb7517SMichael Lange     }
342331830f3SMatthew G. Knepley   }
3436228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
344c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
345331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
346331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
347331830f3SMatthew G. Knepley   if (interpolate) {
348e56d480eSMatthew G. Knepley     DM idm = NULL;
349331830f3SMatthew G. Knepley 
350331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
351331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
352331830f3SMatthew G. Knepley     *dm  = idm;
353331830f3SMatthew G. Knepley   }
35416ad7e67SMichael Lange 
35516ad7e67SMichael Lange   if (!rank) {
35616ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
35716ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
358d1a54cefSMatthew G. Knepley 
35916ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
36016ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
361a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
36216ad7e67SMichael Lange         PetscInt joinSize;
36316ad7e67SMichael Lange         const PetscInt *join;
364a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
365a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
36616ad7e67SMichael Lange         }
367a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
368a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
369c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
370a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
37116ad7e67SMichael Lange       }
372a4bb7517SMichael Lange     }
3736e1dd89aSLawrence Mitchell 
3746e1dd89aSLawrence Mitchell     /* Create cell sets */
3756e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
3766e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
3776e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
3786e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
3796e1dd89aSLawrence Mitchell           cell++;
3806e1dd89aSLawrence Mitchell         }
3816e1dd89aSLawrence Mitchell       }
3826e1dd89aSLawrence Mitchell     }
3830c070f12SSander Arens 
3840c070f12SSander Arens     /* Create vertex sets */
3850c070f12SSander Arens     for (c = 0; c < numCells; ++c) {
3860c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
3870c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
3880c070f12SSander Arens           ierr = DMSetLabelValue(*dm, "Vertex Sets", gmsh_elem[c].nodes[0] + vStart - 1, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
3890c070f12SSander Arens         }
3900c070f12SSander Arens       }
3910c070f12SSander Arens     }
39216ad7e67SMichael Lange   }
39316ad7e67SMichael Lange 
394331830f3SMatthew G. Knepley   /* Read coordinates */
395979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
396331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
397331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
39875b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
39975b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
400331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
401331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
402331830f3SMatthew G. Knepley   }
403331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
404331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
4058b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
406331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
407331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4087635fff5Ssarens   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
409331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
410331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
411331830f3SMatthew G. Knepley   if (!rank) {
412331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
413331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
414331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
415331830f3SMatthew G. Knepley       }
416331830f3SMatthew G. Knepley     }
417331830f3SMatthew G. Knepley   }
418331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
419331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
420331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
421331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
422a4bb7517SMichael Lange   /* Clean up intermediate storage */
42329403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
4243b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
425331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
426331830f3SMatthew G. Knepley }
427