xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 21016a8b907ae59b2cfc4d906ff16c10bbb7e7f3)
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;
91*21016a8bSBarry Smith   PetscScalar     *coords;
92*21016a8bSBarry Smith   PetscReal       *coordsIn = NULL;
93dc0ede3bSMatthew G. Knepley   PetscInt        dim = 0, coordSize, c, v, d, r, cell;
94dc0ede3bSMatthew G. Knepley   int             i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
95331830f3SMatthew G. Knepley   PetscMPIInt     num_proc, rank;
96331830f3SMatthew G. Knepley   char            line[PETSC_MAX_PATH_LEN];
9704d1ad83SMichael Lange   PetscBool       match, binary, bswap = PETSC_FALSE;
98331830f3SMatthew G. Knepley   PetscErrorCode  ierr;
99331830f3SMatthew G. Knepley 
100331830f3SMatthew G. Knepley   PetscFunctionBegin;
101331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
102331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
103331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
104331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1053b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
10604d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
10704d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
108abc86ac4SMichael Lange   if (!rank || binary) {
109331830f3SMatthew G. Knepley     PetscBool match;
110331830f3SMatthew G. Knepley     int       fileType, dataSize;
111f6ac7a6aSMichael Lange     float     version;
112331830f3SMatthew G. Knepley 
113331830f3SMatthew G. Knepley     /* Read header */
114060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
11504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
116331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
117060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
118f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
119f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
120f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
121331830f3SMatthew 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);
12204d1ad83SMichael Lange     if (binary) {
123b9eae255SMichael Lange       int checkInt;
124060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
12504d1ad83SMichael Lange       if (checkInt != 1) {
126b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
12704d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
12804d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
12904d1ad83SMichael Lange       }
1300877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
131060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
13204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
133331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
134dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
135060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
136dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
137dc0ede3bSMatthew G. Knepley     if (match) {
138dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
139dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
140dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
141a8667449SSander Arens       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
142dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
143dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
144dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
145dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
146dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
147dc0ede3bSMatthew G. Knepley     }
148dc0ede3bSMatthew G. Knepley     /* Read vertices */
14904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
150331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
151060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
1520877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
1530877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
154854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
15513aecfe5SMichael Lange     if (binary) {
15613aecfe5SMichael Lange       size_t doubleSize, intSize;
15713aecfe5SMichael Lange       PetscInt elementSize;
15813aecfe5SMichael Lange       char *buffer;
15913aecfe5SMichael Lange       PetscScalar *baseptr;
16013aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
16113aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
16213aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
16313aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
16413aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
16513aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
166331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
16713aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
16813aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
16913aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
17013aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
17113aecfe5SMichael Lange       }
17213aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
17313aecfe5SMichael Lange     } else {
17413aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
175060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
176fba955ccSBarry Smith         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_REAL);CHKERRQ(ierr);
177fba955ccSBarry Smith         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %D should be %D", i, v+1);
178331830f3SMatthew G. Knepley       }
17913aecfe5SMichael Lange     }
180060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
18104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
182331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
183331830f3SMatthew G. Knepley     /* Read cells */
184060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
18504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
186331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
187060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
1880877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
1890877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
190331830f3SMatthew G. Knepley   }
191db397164SMichael Lange 
192abc86ac4SMichael Lange   if (!rank || binary) {
193a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
194a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
195a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
196a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
1973d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
198a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
199a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
200a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
201db397164SMichael Lange     }
202060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2031b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
204331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
205331830f3SMatthew G. Knepley   }
206abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
207abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
208a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
209331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
210331830f3SMatthew G. Knepley   if (!rank) {
211a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
212a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
213a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
214a4bb7517SMichael Lange         cell++;
215331830f3SMatthew G. Knepley       }
216331830f3SMatthew G. Knepley     }
217331830f3SMatthew G. Knepley   }
218331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
219a4bb7517SMichael Lange   /* Add cell-vertex connections */
220331830f3SMatthew G. Knepley   if (!rank) {
221331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
222a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
223a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
224a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
225a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1;
226331830f3SMatthew G. Knepley         }
22797ed6be6Ssarens         if (dim == 3) {
22897ed6be6Ssarens           /* Tetrahedra are inverted */
22997ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
23097ed6be6Ssarens             PetscInt tmp = pcone[0];
23197ed6be6Ssarens             pcone[0] = pcone[1];
23297ed6be6Ssarens             pcone[1] = tmp;
23397ed6be6Ssarens           }
23497ed6be6Ssarens           /* Hexahedra are inverted */
23597ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
23697ed6be6Ssarens             PetscInt tmp = pcone[1];
23797ed6be6Ssarens             pcone[1] = pcone[3];
23897ed6be6Ssarens             pcone[3] = tmp;
23997ed6be6Ssarens           }
24097ed6be6Ssarens         }
241a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
242a4bb7517SMichael Lange         cell++;
243331830f3SMatthew G. Knepley       }
244a4bb7517SMichael Lange     }
245331830f3SMatthew G. Knepley   }
2466228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
247c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
248331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
249331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
250331830f3SMatthew G. Knepley   if (interpolate) {
251e56d480eSMatthew G. Knepley     DM idm = NULL;
252331830f3SMatthew G. Knepley 
253331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
254331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
255331830f3SMatthew G. Knepley     *dm  = idm;
256331830f3SMatthew G. Knepley   }
25716ad7e67SMichael Lange 
25816ad7e67SMichael Lange   if (!rank) {
25916ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
26016ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
261d1a54cefSMatthew G. Knepley 
26216ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
26316ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
264a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
26516ad7e67SMichael Lange         PetscInt joinSize;
26616ad7e67SMichael Lange         const PetscInt *join;
267a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
268a4bb7517SMichael Lange           pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1;
26916ad7e67SMichael Lange         }
270a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
271a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
272c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
273a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
27416ad7e67SMichael Lange       }
275a4bb7517SMichael Lange     }
2766e1dd89aSLawrence Mitchell 
2776e1dd89aSLawrence Mitchell     /* Create cell sets */
2786e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
2796e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
2806e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
2816e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
2826e1dd89aSLawrence Mitchell           cell++;
2836e1dd89aSLawrence Mitchell         }
2846e1dd89aSLawrence Mitchell       }
2856e1dd89aSLawrence Mitchell     }
2860c070f12SSander Arens 
2870c070f12SSander Arens     /* Create vertex sets */
2880c070f12SSander Arens     for (c = 0; c < numCells; ++c) {
2890c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
2900c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
2910c070f12SSander Arens           ierr = DMSetLabelValue(*dm, "Vertex Sets", gmsh_elem[c].nodes[0] + vStart - 1, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
2920c070f12SSander Arens         }
2930c070f12SSander Arens       }
2940c070f12SSander Arens     }
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);
3088b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
309331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
310331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
3117635fff5Ssarens   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 
3313d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
332db397164SMichael Lange {
3335257f852SMichael Lange   PetscInt       c, p;
3343d51f2daSMichael Lange   GmshElement   *elements;
335b6a3fe3cSMichael Lange   int            i, cellType, dim, numNodes, numElem, numTags;
3365257f852SMichael Lange   int            ibuf[16];
337a4bb7517SMichael Lange   PetscErrorCode ierr;
338db397164SMichael Lange 
33930412aabSMichael Lange   PetscFunctionBegin;
3403d51f2daSMichael Lange   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
3413d51f2daSMichael Lange   for (c = 0; c < numCells;) {
342b6a3fe3cSMichael Lange     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
343b6a3fe3cSMichael Lange     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
34404d1ad83SMichael Lange     if (binary) {
345b6a3fe3cSMichael Lange       cellType = ibuf[0];
346b6a3fe3cSMichael Lange       numElem = ibuf[1];
347b6a3fe3cSMichael Lange       numTags = ibuf[2];
34804d1ad83SMichael Lange     } else {
349b6a3fe3cSMichael Lange       elements[c].id = ibuf[0];
350b6a3fe3cSMichael Lange       cellType = ibuf[1];
351b6a3fe3cSMichael Lange       numTags = ibuf[2];
3523d51f2daSMichael Lange       numElem = 1;
35304d1ad83SMichael Lange     }
354b6a3fe3cSMichael Lange     switch (cellType) {
355b6a3fe3cSMichael Lange     case 1: /* 2-node line */
356b6a3fe3cSMichael Lange       dim = 1;
357b6a3fe3cSMichael Lange       numNodes = 2;
358b6a3fe3cSMichael Lange       break;
359b6a3fe3cSMichael Lange     case 2: /* 3-node triangle */
360b6a3fe3cSMichael Lange       dim = 2;
361b6a3fe3cSMichael Lange       numNodes = 3;
362b6a3fe3cSMichael Lange       break;
363b6a3fe3cSMichael Lange     case 3: /* 4-node quadrangle */
364b6a3fe3cSMichael Lange       dim = 2;
365b6a3fe3cSMichael Lange       numNodes = 4;
366b6a3fe3cSMichael Lange       break;
367b6a3fe3cSMichael Lange     case 4: /* 4-node tetrahedron */
368b6a3fe3cSMichael Lange       dim  = 3;
369b6a3fe3cSMichael Lange       numNodes = 4;
370b6a3fe3cSMichael Lange       break;
371b6a3fe3cSMichael Lange     case 5: /* 8-node hexahedron */
372b6a3fe3cSMichael Lange       dim = 3;
373b6a3fe3cSMichael Lange       numNodes = 8;
374b6a3fe3cSMichael Lange       break;
375b6a3fe3cSMichael Lange     case 15: /* 1-node vertex */
376b6a3fe3cSMichael Lange       dim = 0;
377b6a3fe3cSMichael Lange       numNodes = 1;
378b6a3fe3cSMichael Lange       break;
379b6a3fe3cSMichael Lange     default:
380b6a3fe3cSMichael Lange       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
381b6a3fe3cSMichael Lange     }
3825257f852SMichael Lange     if (binary) {
3835257f852SMichael Lange       const PetscInt nint = numNodes + numTags + 1;
3843d51f2daSMichael Lange       for (i = 0; i < numElem; ++i, ++c) {
3853d51f2daSMichael Lange         /* Loop over inner binary element block */
386b6a3fe3cSMichael Lange         elements[c].dim = dim;
387b6a3fe3cSMichael Lange         elements[c].numNodes = numNodes;
388b6a3fe3cSMichael Lange         elements[c].numTags = numTags;
389b6a3fe3cSMichael Lange 
3905257f852SMichael Lange         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
3915257f852SMichael Lange         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
3925257f852SMichael Lange         elements[c].id = ibuf[0];
3935257f852SMichael Lange         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
3945257f852SMichael Lange         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
3955257f852SMichael Lange       }
3965257f852SMichael Lange     } else {
3975257f852SMichael Lange       elements[c].dim = dim;
3985257f852SMichael Lange       elements[c].numNodes = numNodes;
3995257f852SMichael Lange       elements[c].numTags = numTags;
4003d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
4013d51f2daSMichael Lange       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
4025257f852SMichael Lange       c++;
4033d51f2daSMichael Lange     }
4043d51f2daSMichael Lange   }
4053d51f2daSMichael Lange   *gmsh_elems = elements;
40630412aabSMichael Lange   PetscFunctionReturn(0);
407db397164SMichael Lange }
408