xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 1b82c3eb2ec338d496a42e9b957dda5e8108f9c8)
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 
637d282ae0SMichael Lange 
64331830f3SMatthew G. Knepley /*@
657d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
66331830f3SMatthew G. Knepley 
67331830f3SMatthew G. Knepley   Collective on comm
68331830f3SMatthew G. Knepley 
69331830f3SMatthew G. Knepley   Input Parameters:
70331830f3SMatthew G. Knepley + comm  - The MPI communicator
71331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
72331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
73331830f3SMatthew G. Knepley 
74331830f3SMatthew G. Knepley   Output Parameter:
75331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
76331830f3SMatthew G. Knepley 
77331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
783d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
79331830f3SMatthew G. Knepley 
80331830f3SMatthew G. Knepley   Level: beginner
81331830f3SMatthew G. Knepley 
82331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
83331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
84331830f3SMatthew G. Knepley @*/
85331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
86331830f3SMatthew G. Knepley {
8704d1ad83SMichael Lange   PetscViewerType vtype;
88a4bb7517SMichael Lange   GmshElement   *gmsh_elem;
89331830f3SMatthew G. Knepley   PetscSection   coordSection;
90331830f3SMatthew G. Knepley   Vec            coordinates;
91331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
92dc0ede3bSMatthew G. Knepley   PetscInt       dim = 0, coordSize, c, v, d, r, cell;
93dc0ede3bSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
94331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
95331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
9604d1ad83SMichael Lange   PetscBool      match, binary, bswap = PETSC_FALSE;
97331830f3SMatthew G. Knepley   PetscErrorCode ierr;
98331830f3SMatthew G. Knepley 
99331830f3SMatthew G. Knepley   PetscFunctionBegin;
100331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
101331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
102331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
103331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1043b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);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;
110f6ac7a6aSMichael Lange     float     version;
111331830f3SMatthew G. Knepley 
112331830f3SMatthew G. Knepley     /* Read header */
113060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, 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");
116060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
117f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
118f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
119f6ac7a6aSMichael 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;
123060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, 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);
130060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, 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");
133dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
134060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
135dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
136dc0ede3bSMatthew G. Knepley     if (match) {
137dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
138dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
139dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
140a8667449SSander Arens       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
141dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
142dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
143dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
144dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
145dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
146dc0ede3bSMatthew G. Knepley     }
147dc0ede3bSMatthew G. Knepley     /* Read vertices */
14804d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
149331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
150060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
1510877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
1520877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
153854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
15413aecfe5SMichael Lange     if (binary) {
15513aecfe5SMichael Lange       size_t doubleSize, intSize;
15613aecfe5SMichael Lange       PetscInt elementSize;
15713aecfe5SMichael Lange       char *buffer;
15813aecfe5SMichael Lange       PetscScalar *baseptr;
15913aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
16013aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
16113aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
16213aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
16313aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
16413aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
165331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
16613aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
16713aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
16813aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
16913aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
17013aecfe5SMichael Lange       }
17113aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
17213aecfe5SMichael Lange     } else {
17313aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
174060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
175060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
176b9eae255SMichael Lange         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
177331830f3SMatthew G. Knepley       }
17813aecfe5SMichael Lange     }
179060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
18004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
181331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
182331830f3SMatthew G. Knepley     /* Read cells */
183060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
18404d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
185331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
186060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
1870877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
1880877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
189331830f3SMatthew G. Knepley   }
190db397164SMichael Lange 
191abc86ac4SMichael Lange   if (!rank || binary) {
192a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
193a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
194a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
195a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
1963d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
197a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
198a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
199a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
200db397164SMichael Lange     }
201060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
202*1b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
203331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
204331830f3SMatthew G. Knepley   }
205abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
206abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
207a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
208331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
209331830f3SMatthew G. Knepley   if (!rank) {
210a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
211a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
212a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
213a4bb7517SMichael Lange         cell++;
214331830f3SMatthew G. Knepley       }
215331830f3SMatthew G. Knepley     }
216331830f3SMatthew G. Knepley   }
217331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
218a4bb7517SMichael Lange   /* Add cell-vertex connections */
219331830f3SMatthew G. Knepley   if (!rank) {
220331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
221a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
222a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
223a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
224a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
225331830f3SMatthew G. Knepley         }
22697ed6be6Ssarens         if (dim == 3) {
22797ed6be6Ssarens           /* Tetrahedra are inverted */
22897ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
22997ed6be6Ssarens             PetscInt tmp = pcone[0];
23097ed6be6Ssarens             pcone[0] = pcone[1];
23197ed6be6Ssarens             pcone[1] = tmp;
23297ed6be6Ssarens           }
23397ed6be6Ssarens           /* Hexahedra are inverted */
23497ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
23597ed6be6Ssarens             PetscInt tmp = pcone[1];
23697ed6be6Ssarens             pcone[1] = pcone[3];
23797ed6be6Ssarens             pcone[3] = tmp;
23897ed6be6Ssarens           }
23997ed6be6Ssarens         }
240a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
241a4bb7517SMichael Lange         cell++;
242331830f3SMatthew G. Knepley       }
243a4bb7517SMichael Lange     }
244331830f3SMatthew G. Knepley   }
2456228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
246c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
247331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
249331830f3SMatthew G. Knepley   if (interpolate) {
250e56d480eSMatthew G. Knepley     DM idm = NULL;
251331830f3SMatthew G. Knepley 
252331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
253331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
254331830f3SMatthew G. Knepley     *dm  = idm;
255331830f3SMatthew G. Knepley   }
25616ad7e67SMichael Lange 
25716ad7e67SMichael Lange   if (!rank) {
25816ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
25916ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
260d1a54cefSMatthew G. Knepley 
26116ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
26216ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
263a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
26416ad7e67SMichael Lange         PetscInt joinSize;
26516ad7e67SMichael Lange         const PetscInt *join;
266a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
267a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
26816ad7e67SMichael Lange         }
269a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
270a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
271c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
272a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
27316ad7e67SMichael Lange       }
274a4bb7517SMichael Lange     }
2756e1dd89aSLawrence Mitchell 
2766e1dd89aSLawrence Mitchell     /* Create cell sets */
2776e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
2786e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
2796e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
2806e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
2816e1dd89aSLawrence Mitchell           cell++;
2826e1dd89aSLawrence Mitchell         }
2836e1dd89aSLawrence Mitchell       }
2846e1dd89aSLawrence Mitchell     }
2850c070f12SSander Arens 
2860c070f12SSander Arens     /* Create vertex sets */
2870c070f12SSander Arens     for (c = 0; c < numCells; ++c) {
2880c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
2890c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
2900c070f12SSander Arens           ierr = DMSetLabelValue(*dm, "Vertex Sets", gmsh_elem[c].nodes[0] + vStart - 1, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
2910c070f12SSander Arens         }
2920c070f12SSander Arens       }
2930c070f12SSander Arens     }
29416ad7e67SMichael Lange   }
29516ad7e67SMichael Lange 
296331830f3SMatthew G. Knepley   /* Read coordinates */
297979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
298331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
299331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
30075b5763bSMichael Lange   ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
30175b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
302331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
303331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
304331830f3SMatthew G. Knepley   }
305331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
306331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
3078b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
308331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
309331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
3107635fff5Ssarens   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
311331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
312331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
313331830f3SMatthew G. Knepley   if (!rank) {
314331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
315331830f3SMatthew G. Knepley       for (d = 0; d < dim; ++d) {
316331830f3SMatthew G. Knepley         coords[v*dim+d] = coordsIn[v*3+d];
317331830f3SMatthew G. Knepley       }
318331830f3SMatthew G. Knepley     }
319331830f3SMatthew G. Knepley   }
320331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
321331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
322331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
323331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
324a4bb7517SMichael Lange   /* Clean up intermediate storage */
32529403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
3263b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
327331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
328331830f3SMatthew G. Knepley }
329db397164SMichael Lange 
3303d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
331db397164SMichael Lange {
3325257f852SMichael Lange   PetscInt       c, p;
3333d51f2daSMichael Lange   GmshElement   *elements;
334b6a3fe3cSMichael Lange   int            i, cellType, dim, numNodes, numElem, numTags;
3355257f852SMichael Lange   int            ibuf[16];
336a4bb7517SMichael Lange   PetscErrorCode ierr;
337db397164SMichael Lange 
33830412aabSMichael Lange   PetscFunctionBegin;
3393d51f2daSMichael Lange   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
3403d51f2daSMichael Lange   for (c = 0; c < numCells;) {
341b6a3fe3cSMichael Lange     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
342b6a3fe3cSMichael Lange     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
34304d1ad83SMichael Lange     if (binary) {
344b6a3fe3cSMichael Lange       cellType = ibuf[0];
345b6a3fe3cSMichael Lange       numElem = ibuf[1];
346b6a3fe3cSMichael Lange       numTags = ibuf[2];
34704d1ad83SMichael Lange     } else {
348b6a3fe3cSMichael Lange       elements[c].id = ibuf[0];
349b6a3fe3cSMichael Lange       cellType = ibuf[1];
350b6a3fe3cSMichael Lange       numTags = ibuf[2];
3513d51f2daSMichael Lange       numElem = 1;
35204d1ad83SMichael Lange     }
353b6a3fe3cSMichael Lange     switch (cellType) {
354b6a3fe3cSMichael Lange     case 1: /* 2-node line */
355b6a3fe3cSMichael Lange       dim = 1;
356b6a3fe3cSMichael Lange       numNodes = 2;
357b6a3fe3cSMichael Lange       break;
358b6a3fe3cSMichael Lange     case 2: /* 3-node triangle */
359b6a3fe3cSMichael Lange       dim = 2;
360b6a3fe3cSMichael Lange       numNodes = 3;
361b6a3fe3cSMichael Lange       break;
362b6a3fe3cSMichael Lange     case 3: /* 4-node quadrangle */
363b6a3fe3cSMichael Lange       dim = 2;
364b6a3fe3cSMichael Lange       numNodes = 4;
365b6a3fe3cSMichael Lange       break;
366b6a3fe3cSMichael Lange     case 4: /* 4-node tetrahedron */
367b6a3fe3cSMichael Lange       dim  = 3;
368b6a3fe3cSMichael Lange       numNodes = 4;
369b6a3fe3cSMichael Lange       break;
370b6a3fe3cSMichael Lange     case 5: /* 8-node hexahedron */
371b6a3fe3cSMichael Lange       dim = 3;
372b6a3fe3cSMichael Lange       numNodes = 8;
373b6a3fe3cSMichael Lange       break;
374b6a3fe3cSMichael Lange     case 15: /* 1-node vertex */
375b6a3fe3cSMichael Lange       dim = 0;
376b6a3fe3cSMichael Lange       numNodes = 1;
377b6a3fe3cSMichael Lange       break;
378b6a3fe3cSMichael Lange     default:
379b6a3fe3cSMichael Lange       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
380b6a3fe3cSMichael Lange     }
3815257f852SMichael Lange     if (binary) {
3825257f852SMichael Lange       const PetscInt nint = numNodes + numTags + 1;
3833d51f2daSMichael Lange       for (i = 0; i < numElem; ++i, ++c) {
3843d51f2daSMichael Lange         /* Loop over inner binary element block */
385b6a3fe3cSMichael Lange         elements[c].dim = dim;
386b6a3fe3cSMichael Lange         elements[c].numNodes = numNodes;
387b6a3fe3cSMichael Lange         elements[c].numTags = numTags;
388b6a3fe3cSMichael Lange 
3895257f852SMichael Lange         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
3905257f852SMichael Lange         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
3915257f852SMichael Lange         elements[c].id = ibuf[0];
3925257f852SMichael Lange         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
3935257f852SMichael Lange         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
3945257f852SMichael Lange       }
3955257f852SMichael Lange     } else {
3965257f852SMichael Lange       elements[c].dim = dim;
3975257f852SMichael Lange       elements[c].numNodes = numNodes;
3985257f852SMichael Lange       elements[c].numTags = numTags;
3993d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
4003d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
4015257f852SMichael Lange       c++;
4023d51f2daSMichael Lange     }
4033d51f2daSMichael Lange   }
4043d51f2daSMichael Lange   *gmsh_elems = elements;
40530412aabSMichael Lange   PetscFunctionReturn(0);
406db397164SMichael Lange }
407