xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 7635fff5e9b6517fa671206ca03092424fdf1493)
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;
96dc0ede3bSMatthew G. Knepley   PetscInt       dim = 0, coordSize, c, v, d, r, cell;
97dc0ede3bSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 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");
137dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
138060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
139dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
140dc0ede3bSMatthew G. Knepley     if (match) {
141dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
142dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
143dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
144dc0ede3bSMatthew G. Knepley       for (r = 0; r < numRegions; ++r) {
145dc0ede3bSMatthew G. Knepley         PetscInt rdim, tag;
146dc0ede3bSMatthew G. Knepley 
147dc0ede3bSMatthew G. Knepley         ierr = PetscViewerRead(viewer, &rdim, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
148dc0ede3bSMatthew G. Knepley         ierr = PetscViewerRead(viewer, &tag,  1, NULL, PETSC_ENUM);CHKERRQ(ierr);
149dc0ede3bSMatthew G. Knepley         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
150dc0ede3bSMatthew G. Knepley       }
151dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
152dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
153dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
154dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
155dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
156dc0ede3bSMatthew G. Knepley     }
157dc0ede3bSMatthew G. Knepley     /* Read vertices */
15804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", 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", &numVertices);
1620877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
163854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
16413aecfe5SMichael Lange     if (binary) {
16513aecfe5SMichael Lange       size_t doubleSize, intSize;
16613aecfe5SMichael Lange       PetscInt elementSize;
16713aecfe5SMichael Lange       char *buffer;
16813aecfe5SMichael Lange       PetscScalar *baseptr;
16913aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
17013aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
17113aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
17213aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
17313aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
17413aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
175331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
17613aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
17713aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
17813aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
17913aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
18013aecfe5SMichael Lange       }
18113aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
18213aecfe5SMichael Lange     } else {
18313aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
184060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
185060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
186b9eae255SMichael Lange         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
187331830f3SMatthew G. Knepley       }
18813aecfe5SMichael Lange     }
189060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
19004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
191331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
192331830f3SMatthew G. Knepley     /* Read cells */
193060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
19404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
195331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
196060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
1970877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
1980877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
199331830f3SMatthew G. Knepley   }
200db397164SMichael Lange 
201abc86ac4SMichael Lange   if (!rank || binary) {
202a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
203a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
204a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
205a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
2063d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
207a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
208a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
209a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
210db397164SMichael Lange     }
211060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
21204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndElements", 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");
214331830f3SMatthew G. Knepley   }
215abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
216abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
217a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
218331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
219331830f3SMatthew G. Knepley   if (!rank) {
220a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
221a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
222a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
223a4bb7517SMichael Lange         cell++;
224331830f3SMatthew G. Knepley       }
225331830f3SMatthew G. Knepley     }
226331830f3SMatthew G. Knepley   }
227331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
228a4bb7517SMichael Lange   /* Add cell-vertex connections */
229331830f3SMatthew G. Knepley   if (!rank) {
230331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
231a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
232a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
233a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
234a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
235331830f3SMatthew G. Knepley         }
23697ed6be6Ssarens         if (dim == 3) {
23797ed6be6Ssarens           /* Tetrahedra are inverted */
23897ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
23997ed6be6Ssarens             PetscInt tmp = pcone[0];
24097ed6be6Ssarens             pcone[0] = pcone[1];
24197ed6be6Ssarens             pcone[1] = tmp;
24297ed6be6Ssarens           }
24397ed6be6Ssarens           /* Hexahedra are inverted */
24497ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
24597ed6be6Ssarens             PetscInt tmp = pcone[1];
24697ed6be6Ssarens             pcone[1] = pcone[3];
24797ed6be6Ssarens             pcone[3] = tmp;
24897ed6be6Ssarens           }
24997ed6be6Ssarens         }
250a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
251a4bb7517SMichael Lange         cell++;
252331830f3SMatthew G. Knepley       }
253a4bb7517SMichael Lange     }
254331830f3SMatthew G. Knepley   }
2556228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
256c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
257331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
258331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
259331830f3SMatthew G. Knepley   if (interpolate) {
260e56d480eSMatthew G. Knepley     DM idm = NULL;
261331830f3SMatthew G. Knepley 
262331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
263331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
264331830f3SMatthew G. Knepley     *dm  = idm;
265331830f3SMatthew G. Knepley   }
26616ad7e67SMichael Lange 
26716ad7e67SMichael Lange   if (!rank) {
26816ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
26916ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
270d1a54cefSMatthew G. Knepley 
27116ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
27216ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
273a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
27416ad7e67SMichael Lange         PetscInt joinSize;
27516ad7e67SMichael Lange         const PetscInt *join;
276a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
277a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
27816ad7e67SMichael Lange         }
279a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
280a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
281c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
282a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
28316ad7e67SMichael Lange       }
284a4bb7517SMichael Lange     }
2856e1dd89aSLawrence Mitchell 
2866e1dd89aSLawrence Mitchell     /* Create cell sets */
2876e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
2886e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
2896e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
2906e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
2916e1dd89aSLawrence Mitchell           cell++;
2926e1dd89aSLawrence Mitchell         }
2936e1dd89aSLawrence Mitchell       }
2946e1dd89aSLawrence Mitchell     }
29516ad7e67SMichael Lange   }
29616ad7e67SMichael Lange 
297331830f3SMatthew G. Knepley   /* Read coordinates */
298979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
299331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
300331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
30175b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
30275b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
303331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
304331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
305331830f3SMatthew G. Knepley   }
306331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
307331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
308331830f3SMatthew G. Knepley   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
309331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
310331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
311*7635fff5Ssarens   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
312331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
313331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
314331830f3SMatthew G. Knepley   if (!rank) {
315331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
316331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
317331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
318331830f3SMatthew G. Knepley       }
319331830f3SMatthew G. Knepley     }
320331830f3SMatthew G. Knepley   }
321331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
322331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
323331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
324331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
325a4bb7517SMichael Lange   /* Clean up intermediate storage */
32629403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
3273b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
328331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
329331830f3SMatthew G. Knepley }
330db397164SMichael Lange 
331db397164SMichael Lange #undef __FUNCT__
33230412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement"
3333d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
334db397164SMichael Lange {
3355257f852SMichael Lange   PetscInt       c, p;
3363d51f2daSMichael Lange   GmshElement   *elements;
337b6a3fe3cSMichael Lange   int            i, cellType, dim, numNodes, numElem, numTags;
3385257f852SMichael Lange   int            ibuf[16];
339a4bb7517SMichael Lange   PetscErrorCode ierr;
340db397164SMichael Lange 
34130412aabSMichael Lange   PetscFunctionBegin;
3423d51f2daSMichael Lange   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
3433d51f2daSMichael Lange   for (c = 0; c < numCells;) {
344b6a3fe3cSMichael Lange     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
345b6a3fe3cSMichael Lange     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
34604d1ad83SMichael Lange     if (binary) {
347b6a3fe3cSMichael Lange       cellType = ibuf[0];
348b6a3fe3cSMichael Lange       numElem = ibuf[1];
349b6a3fe3cSMichael Lange       numTags = ibuf[2];
35004d1ad83SMichael Lange     } else {
351b6a3fe3cSMichael Lange       elements[c].id = ibuf[0];
352b6a3fe3cSMichael Lange       cellType = ibuf[1];
353b6a3fe3cSMichael Lange       numTags = ibuf[2];
3543d51f2daSMichael Lange       numElem = 1;
35504d1ad83SMichael Lange     }
356b6a3fe3cSMichael Lange     switch (cellType) {
357b6a3fe3cSMichael Lange     case 1: /* 2-node line */
358b6a3fe3cSMichael Lange       dim = 1;
359b6a3fe3cSMichael Lange       numNodes = 2;
360b6a3fe3cSMichael Lange       break;
361b6a3fe3cSMichael Lange     case 2: /* 3-node triangle */
362b6a3fe3cSMichael Lange       dim = 2;
363b6a3fe3cSMichael Lange       numNodes = 3;
364b6a3fe3cSMichael Lange       break;
365b6a3fe3cSMichael Lange     case 3: /* 4-node quadrangle */
366b6a3fe3cSMichael Lange       dim = 2;
367b6a3fe3cSMichael Lange       numNodes = 4;
368b6a3fe3cSMichael Lange       break;
369b6a3fe3cSMichael Lange     case 4: /* 4-node tetrahedron */
370b6a3fe3cSMichael Lange       dim  = 3;
371b6a3fe3cSMichael Lange       numNodes = 4;
372b6a3fe3cSMichael Lange       break;
373b6a3fe3cSMichael Lange     case 5: /* 8-node hexahedron */
374b6a3fe3cSMichael Lange       dim = 3;
375b6a3fe3cSMichael Lange       numNodes = 8;
376b6a3fe3cSMichael Lange       break;
377b6a3fe3cSMichael Lange     case 15: /* 1-node vertex */
378b6a3fe3cSMichael Lange       dim = 0;
379b6a3fe3cSMichael Lange       numNodes = 1;
380b6a3fe3cSMichael Lange       break;
381b6a3fe3cSMichael Lange     default:
382b6a3fe3cSMichael Lange       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
383b6a3fe3cSMichael Lange     }
3845257f852SMichael Lange     if (binary) {
3855257f852SMichael Lange       const PetscInt nint = numNodes + numTags + 1;
3863d51f2daSMichael Lange       for (i = 0; i < numElem; ++i, ++c) {
3873d51f2daSMichael Lange         /* Loop over inner binary element block */
388b6a3fe3cSMichael Lange         elements[c].dim = dim;
389b6a3fe3cSMichael Lange         elements[c].numNodes = numNodes;
390b6a3fe3cSMichael Lange         elements[c].numTags = numTags;
391b6a3fe3cSMichael Lange 
3925257f852SMichael Lange         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
3935257f852SMichael Lange         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
3945257f852SMichael Lange         elements[c].id = ibuf[0];
3955257f852SMichael Lange         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
3965257f852SMichael Lange         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
3975257f852SMichael Lange       }
3985257f852SMichael Lange     } else {
3995257f852SMichael Lange       elements[c].dim = dim;
4005257f852SMichael Lange       elements[c].numNodes = numNodes;
4015257f852SMichael Lange       elements[c].numTags = numTags;
4023d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
4033d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
4045257f852SMichael Lange       c++;
4053d51f2daSMichael Lange     }
4063d51f2daSMichael Lange   }
4073d51f2daSMichael Lange   *gmsh_elems = elements;
40830412aabSMichael Lange   PetscFunctionReturn(0);
409db397164SMichael Lange }
410