xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 730356f691fef1270aa21ce63f6c07feb0212308)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3*730356f6SLisandro Dalcin #include <petsc/private/hashmapi.h>
4331830f3SMatthew G. Knepley 
57d282ae0SMichael Lange /*@C
67d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
77d282ae0SMichael Lange 
87d282ae0SMichael Lange + comm        - The MPI communicator
97d282ae0SMichael Lange . filename    - Name of the Gmsh file
107d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
117d282ae0SMichael Lange 
127d282ae0SMichael Lange   Output Parameter:
137d282ae0SMichael Lange . dm  - The DM object representing the mesh
147d282ae0SMichael Lange 
157d282ae0SMichael Lange   Level: beginner
167d282ae0SMichael Lange 
177d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
187d282ae0SMichael Lange @*/
197d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
207d282ae0SMichael Lange {
2111c56cb0SLisandro Dalcin   PetscViewer     viewer;
22abc86ac4SMichael Lange   PetscMPIInt     rank;
2311c56cb0SLisandro Dalcin   int             fileType;
247d282ae0SMichael Lange   PetscViewerType vtype;
257d282ae0SMichael Lange   PetscErrorCode  ierr;
267d282ae0SMichael Lange 
277d282ae0SMichael Lange   PetscFunctionBegin;
28abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2911c56cb0SLisandro Dalcin 
307d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
3111c56cb0SLisandro Dalcin   if (!rank) {
3211c56cb0SLisandro Dalcin     PetscViewer vheader;
3311c56cb0SLisandro Dalcin     char        line[PETSC_MAX_PATH_LEN];
3411c56cb0SLisandro Dalcin     PetscBool   match;
3511c56cb0SLisandro Dalcin     int         snum;
3611c56cb0SLisandro Dalcin     float       version;
3711c56cb0SLisandro Dalcin 
3811c56cb0SLisandro Dalcin     ierr = PetscViewerCreate(PETSC_COMM_SELF, &vheader);CHKERRQ(ierr);
397d282ae0SMichael Lange     ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
407d282ae0SMichael Lange     ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
417d282ae0SMichael Lange     ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
427d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
43060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
4411c56cb0SLisandro Dalcin     ierr = PetscStrncmp(line, "$MeshFormat", sizeof(line), &match);CHKERRQ(ierr);
457d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
46060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
4711c56cb0SLisandro Dalcin     snum = sscanf(line, "%f %d", &version, &fileType);
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");
50de78e4feSLisandro Dalcin     if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported");
51de78e4feSLisandro Dalcin     if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0");
5211c56cb0SLisandro Dalcin     ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
53abc86ac4SMichael Lange   }
5411c56cb0SLisandro Dalcin   ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
5511c56cb0SLisandro Dalcin   vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
5611c56cb0SLisandro Dalcin 
577d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
587d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
607d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
617d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
627d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
637d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
647d282ae0SMichael Lange   PetscFunctionReturn(0);
657d282ae0SMichael Lange }
667d282ae0SMichael Lange 
67de78e4feSLisandro Dalcin typedef struct {
68de78e4feSLisandro Dalcin   size_t size;
69de78e4feSLisandro Dalcin   void   *buffer;
70de78e4feSLisandro Dalcin } GmshWorkBuffer;
71de78e4feSLisandro Dalcin 
72de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferInit(GmshWorkBuffer *ctx)
73df032b82SMatthew G. Knepley {
74de78e4feSLisandro Dalcin   PetscFunctionBegin;
75de78e4feSLisandro Dalcin   ctx->size   = 0;
76de78e4feSLisandro Dalcin   ctx->buffer = NULL;
77de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
78de78e4feSLisandro Dalcin }
79de78e4feSLisandro Dalcin 
80de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferGet(GmshWorkBuffer *ctx, size_t count, size_t eltsize, void *buf)
81de78e4feSLisandro Dalcin {
82de78e4feSLisandro Dalcin   size_t         size = count*eltsize;
8311c56cb0SLisandro Dalcin   PetscErrorCode ierr;
8411c56cb0SLisandro Dalcin 
8511c56cb0SLisandro Dalcin   PetscFunctionBegin;
86de78e4feSLisandro Dalcin   if (ctx->size < size) {
87de78e4feSLisandro Dalcin     ierr = PetscFree(ctx->buffer);CHKERRQ(ierr);
88de78e4feSLisandro Dalcin     ierr = PetscMalloc(size, &ctx->buffer);CHKERRQ(ierr);
89de78e4feSLisandro Dalcin     ctx->size = size;
90de78e4feSLisandro Dalcin   }
91de78e4feSLisandro Dalcin   *(void**)buf = size ? ctx->buffer : NULL;
92de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
93de78e4feSLisandro Dalcin }
94de78e4feSLisandro Dalcin 
95de78e4feSLisandro Dalcin static PetscErrorCode GmshWorkBufferFree(GmshWorkBuffer *ctx)
96de78e4feSLisandro Dalcin {
97de78e4feSLisandro Dalcin   PetscErrorCode ierr;
98de78e4feSLisandro Dalcin   PetscFunctionBegin;
99de78e4feSLisandro Dalcin   ierr = PetscFree(ctx->buffer);CHKERRQ(ierr);
100de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
101de78e4feSLisandro Dalcin }
102de78e4feSLisandro Dalcin 
103de78e4feSLisandro Dalcin typedef struct {
104de78e4feSLisandro Dalcin   PetscInt dim;      /* Entity dimension */
105de78e4feSLisandro Dalcin   PetscInt id;       /* Entity number */
106de78e4feSLisandro Dalcin   double   bbox[6];  /* Bounding box */
107de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
108de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
109de78e4feSLisandro Dalcin } GmshEntity;
110de78e4feSLisandro Dalcin 
111de78e4feSLisandro Dalcin typedef struct {
112*730356f6SLisandro Dalcin   GmshEntity *entity[4];
113*730356f6SLisandro Dalcin   PetscHMapI  entityMap[4];
114*730356f6SLisandro Dalcin } GmshEntities;
115*730356f6SLisandro Dalcin 
116*730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(long count[4], GmshEntities **entities)
117*730356f6SLisandro Dalcin {
118*730356f6SLisandro Dalcin   PetscInt       dim;
119*730356f6SLisandro Dalcin   PetscErrorCode ierr;
120*730356f6SLisandro Dalcin 
121*730356f6SLisandro Dalcin   PetscFunctionBegin;
122*730356f6SLisandro Dalcin   ierr = PetscNew(entities);CHKERRQ(ierr);
123*730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
124*730356f6SLisandro Dalcin     ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr);
125*730356f6SLisandro Dalcin     ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
126*730356f6SLisandro Dalcin   }
127*730356f6SLisandro Dalcin   PetscFunctionReturn(0);
128*730356f6SLisandro Dalcin }
129*730356f6SLisandro Dalcin 
130*730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity)
131*730356f6SLisandro Dalcin {
132*730356f6SLisandro Dalcin   PetscErrorCode ierr;
133*730356f6SLisandro Dalcin   PetscFunctionBegin;
134*730356f6SLisandro Dalcin   ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr);
135*730356f6SLisandro Dalcin   entities->entity[dim][index].dim = dim;
136*730356f6SLisandro Dalcin   entities->entity[dim][index].id  = eid;
137*730356f6SLisandro Dalcin   if (entity) *entity = &entities->entity[dim][index];
138*730356f6SLisandro Dalcin   PetscFunctionReturn(0);
139*730356f6SLisandro Dalcin }
140*730356f6SLisandro Dalcin 
141*730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity)
142*730356f6SLisandro Dalcin {
143*730356f6SLisandro Dalcin   PetscInt       index;
144*730356f6SLisandro Dalcin   PetscErrorCode ierr;
145*730356f6SLisandro Dalcin 
146*730356f6SLisandro Dalcin   PetscFunctionBegin;
147*730356f6SLisandro Dalcin   ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr);
148*730356f6SLisandro Dalcin   *entity = &entities->entity[dim][index];
149*730356f6SLisandro Dalcin   PetscFunctionReturn(0);
150*730356f6SLisandro Dalcin }
151*730356f6SLisandro Dalcin 
152*730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities)
153*730356f6SLisandro Dalcin {
154*730356f6SLisandro Dalcin   PetscInt       dim;
155*730356f6SLisandro Dalcin   PetscErrorCode ierr;
156*730356f6SLisandro Dalcin 
157*730356f6SLisandro Dalcin   PetscFunctionBegin;
158*730356f6SLisandro Dalcin   if (!*entities) PetscFunctionReturn(0);
159*730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
160*730356f6SLisandro Dalcin     ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr);
161*730356f6SLisandro Dalcin     ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
162*730356f6SLisandro Dalcin   }
163*730356f6SLisandro Dalcin   ierr = PetscFree((*entities));CHKERRQ(ierr);
164*730356f6SLisandro Dalcin   PetscFunctionReturn(0);
165*730356f6SLisandro Dalcin }
166*730356f6SLisandro Dalcin 
167*730356f6SLisandro Dalcin typedef struct {
168de78e4feSLisandro Dalcin   PetscInt dim;      /* Entity dimension */
169de78e4feSLisandro Dalcin   PetscInt id;       /* Element number */
170de78e4feSLisandro Dalcin   PetscInt cellType; /* Cell type */
171de78e4feSLisandro Dalcin   PetscInt numNodes; /* Size of node array */
172de78e4feSLisandro Dalcin   int      nodes[8]; /* Node array */
173de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
174de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
175de78e4feSLisandro Dalcin } GmshElement;
176de78e4feSLisandro Dalcin 
177de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **coordinates)
178de78e4feSLisandro Dalcin {
179de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
180de78e4feSLisandro Dalcin   int            v, nid, snum;
181de78e4feSLisandro Dalcin   PetscErrorCode ierr;
182de78e4feSLisandro Dalcin 
183de78e4feSLisandro Dalcin   PetscFunctionBegin;
184de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
185de78e4feSLisandro Dalcin   snum = sscanf(line, "%d", numVertices);
186de78e4feSLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
187de78e4feSLisandro Dalcin   ierr = PetscMalloc1(*numVertices*3, coordinates);CHKERRQ(ierr);
188de78e4feSLisandro Dalcin   for (v = 0; v < *numVertices; ++v) {
18911c56cb0SLisandro Dalcin     double *xyz = *coordinates + v*3;
19011c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
19111c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
19211c56cb0SLisandro Dalcin     if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift);
19311c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
19411c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
19511c56cb0SLisandro Dalcin   }
19611c56cb0SLisandro Dalcin   PetscFunctionReturn(0);
19711c56cb0SLisandro Dalcin }
19811c56cb0SLisandro Dalcin 
199de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
200de78e4feSLisandro Dalcin    file contents multiple times to figure out the true number of cells and facets
201de78e4feSLisandro Dalcin    in the given mesh. To make this more efficient we read the file contents only
202de78e4feSLisandro Dalcin    once and store them in memory, while determining the true number of cells. */
203de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int *numCells, GmshElement **gmsh_elems)
20411c56cb0SLisandro Dalcin {
205de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
206df032b82SMatthew G. Knepley   GmshElement   *elements;
207de78e4feSLisandro Dalcin   int            i, c, p, ibuf[1+4+512], snum;
20811c56cb0SLisandro Dalcin   int            cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
209df032b82SMatthew G. Knepley   PetscErrorCode ierr;
210df032b82SMatthew G. Knepley 
211df032b82SMatthew G. Knepley   PetscFunctionBegin;
212de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
213de78e4feSLisandro Dalcin   snum = sscanf(line, "%d", numCells);
214de78e4feSLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
215de78e4feSLisandro Dalcin   ierr = PetscMalloc1(*numCells, &elements);CHKERRQ(ierr);
216de78e4feSLisandro Dalcin   for (c = 0; c < *numCells;) {
21711c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
21811c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
219df032b82SMatthew G. Knepley     if (binary) {
220df032b82SMatthew G. Knepley       cellType = ibuf[0];
221df032b82SMatthew G. Knepley       numElem = ibuf[1];
222df032b82SMatthew G. Knepley       numTags = ibuf[2];
223df032b82SMatthew G. Knepley     } else {
224df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
225df032b82SMatthew G. Knepley       cellType = ibuf[1];
226df032b82SMatthew G. Knepley       numTags = ibuf[2];
227df032b82SMatthew G. Knepley       numElem = 1;
228df032b82SMatthew G. Knepley     }
229bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
230bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
231df032b82SMatthew G. Knepley     switch (cellType) {
232df032b82SMatthew G. Knepley     case 1: /* 2-node line */
233df032b82SMatthew G. Knepley       dim = 1;
234df032b82SMatthew G. Knepley       numNodes = 2;
235df032b82SMatthew G. Knepley       break;
236df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
237df032b82SMatthew G. Knepley       dim = 2;
238df032b82SMatthew G. Knepley       numNodes = 3;
239df032b82SMatthew G. Knepley       break;
240df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
241df032b82SMatthew G. Knepley       dim = 2;
242df032b82SMatthew G. Knepley       numNodes = 4;
243df032b82SMatthew G. Knepley       break;
244df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
245df032b82SMatthew G. Knepley       dim  = 3;
246df032b82SMatthew G. Knepley       numNodes = 4;
247df032b82SMatthew G. Knepley       break;
248df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
249df032b82SMatthew G. Knepley       dim = 3;
250df032b82SMatthew G. Knepley       numNodes = 8;
251df032b82SMatthew G. Knepley       break;
252dea2a3c7SStefano Zampini     case 6: /* 6-node wedge */
253dea2a3c7SStefano Zampini       dim = 3;
254dea2a3c7SStefano Zampini       numNodes = 6;
255dea2a3c7SStefano Zampini       break;
256bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
257bf6ba3a3SMatthew G. Knepley       dim = 1;
258bf6ba3a3SMatthew G. Knepley       numNodes = 2;
259bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
260bf6ba3a3SMatthew G. Knepley       break;
261bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
262bf6ba3a3SMatthew G. Knepley       dim = 2;
263bf6ba3a3SMatthew G. Knepley       numNodes = 3;
264bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
265bf6ba3a3SMatthew G. Knepley       break;
266dea2a3c7SStefano Zampini     case 13: /* 18-node 2nd wedge */
267dea2a3c7SStefano Zampini       dim = 3;
268dea2a3c7SStefano Zampini       numNodes = 6;
269dea2a3c7SStefano Zampini       numNodesIgnore = 12;
270dea2a3c7SStefano Zampini       break;
271df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
272df032b82SMatthew G. Knepley       dim = 0;
273df032b82SMatthew G. Knepley       numNodes = 1;
274df032b82SMatthew G. Knepley       break;
275bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
276bf6ba3a3SMatthew G. Knepley     case 10: /* 9-node 2nd order quadrangle */
277bf6ba3a3SMatthew G. Knepley     case 11: /* 10-node 2nd order tetrahedron */
278bf6ba3a3SMatthew G. Knepley     case 12: /* 27-node 2nd order hexhedron */
279bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
280df032b82SMatthew G. Knepley     default:
281df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
282df032b82SMatthew G. Knepley     }
283df032b82SMatthew G. Knepley     if (binary) {
28411c56cb0SLisandro Dalcin       const int nint = 1 + numTags + numNodes + numNodesIgnore;
28511c56cb0SLisandro Dalcin       /* Loop over element blocks */
286df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
28711c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
28811c56cb0SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);}
289df032b82SMatthew G. Knepley         elements[c].dim = dim;
290df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
291df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
292df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
293dea2a3c7SStefano Zampini         elements[c].cellType = cellType;
294df032b82SMatthew G. Knepley         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
295df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
296df032b82SMatthew G. Knepley       }
297df032b82SMatthew G. Knepley     } else {
298df032b82SMatthew G. Knepley       elements[c].dim = dim;
299df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
300df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
301dea2a3c7SStefano Zampini       elements[c].cellType = cellType;
302df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
303df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
30411c56cb0SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
305df032b82SMatthew G. Knepley       c++;
306df032b82SMatthew G. Knepley     }
307df032b82SMatthew G. Knepley   }
308df032b82SMatthew G. Knepley   *gmsh_elems = elements;
309df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
310df032b82SMatthew G. Knepley }
3117d282ae0SMichael Lange 
312de78e4feSLisandro Dalcin /*
313de78e4feSLisandro Dalcin $Entities
314de78e4feSLisandro Dalcin   numPoints(unsigned long) numCurves(unsigned long)
315de78e4feSLisandro Dalcin     numSurfaces(unsigned long) numVolumes(unsigned long)
316de78e4feSLisandro Dalcin   // points
317de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
318de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
319de78e4feSLisandro Dalcin     numPhysicals(unsigned long) phyisicalTag[...](int)
320de78e4feSLisandro Dalcin   ...
321de78e4feSLisandro Dalcin   // curves
322de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
323de78e4feSLisandro Dalcin      boxMaxX(double) boxMaxY(double) boxMaxZ(double)
324de78e4feSLisandro Dalcin      numPhysicals(unsigned long) physicalTag[...](int)
325de78e4feSLisandro Dalcin      numBREPVert(unsigned long) tagBREPVert[...](int)
326de78e4feSLisandro Dalcin   ...
327de78e4feSLisandro Dalcin   // surfaces
328de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
329de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
330de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
331de78e4feSLisandro Dalcin     numBREPCurve(unsigned long) tagBREPCurve[...](int)
332de78e4feSLisandro Dalcin   ...
333de78e4feSLisandro Dalcin   // volumes
334de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
335de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
336de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
337de78e4feSLisandro Dalcin     numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int)
338de78e4feSLisandro Dalcin   ...
339de78e4feSLisandro Dalcin $EndEntities
340de78e4feSLisandro Dalcin */
341*730356f6SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshEntities **entities)
342de78e4feSLisandro Dalcin {
343*730356f6SLisandro Dalcin   long           index, num, count[4];
344*730356f6SLisandro Dalcin   int            dim, eid, numTags, *ibuf, t;
345*730356f6SLisandro Dalcin   GmshEntity     *entity;
346de78e4feSLisandro Dalcin   GmshWorkBuffer work;
347de78e4feSLisandro Dalcin   PetscErrorCode ierr;
348de78e4feSLisandro Dalcin 
349de78e4feSLisandro Dalcin   PetscFunctionBegin;
350de78e4feSLisandro Dalcin   ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr);
351*730356f6SLisandro Dalcin   ierr = PetscViewerRead(viewer, count, 4, NULL, PETSC_LONG);CHKERRQ(ierr);
352*730356f6SLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(count, PETSC_LONG, 4);CHKERRQ(ierr);}
353*730356f6SLisandro Dalcin   ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr);
354de78e4feSLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
355*730356f6SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
356*730356f6SLisandro Dalcin       ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
357*730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);}
358*730356f6SLisandro Dalcin       ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr);
359de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
360de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);}
361de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
362de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
363de78e4feSLisandro Dalcin       ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr);
364de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
365*730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
366de78e4feSLisandro Dalcin       entity->numTags = numTags = (int) PetscMin(num, 4);
367de78e4feSLisandro Dalcin       for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t];
368de78e4feSLisandro Dalcin       if (dim == 0) continue;
369de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
370de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
371de78e4feSLisandro Dalcin       ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr);
372de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
373*730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
374de78e4feSLisandro Dalcin     }
375de78e4feSLisandro Dalcin   }
376de78e4feSLisandro Dalcin   ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr);
377de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
378de78e4feSLisandro Dalcin }
379de78e4feSLisandro Dalcin 
380de78e4feSLisandro Dalcin /*
381de78e4feSLisandro Dalcin $Nodes
382de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numNodes(unsigned long)
383de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long)
384de78e4feSLisandro Dalcin     tag(int) x(double) y(double) z(double)
385de78e4feSLisandro Dalcin     ...
386de78e4feSLisandro Dalcin   ...
387de78e4feSLisandro Dalcin $EndNodes
388de78e4feSLisandro Dalcin */
389de78e4feSLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **gmsh_nodes)
390de78e4feSLisandro Dalcin {
391de78e4feSLisandro Dalcin   long           block, node, v, numEntityBlocks, numTotalNodes, numNodes;
392de78e4feSLisandro Dalcin   int            info[3], nid;
393de78e4feSLisandro Dalcin   double         *coordinates;
394de78e4feSLisandro Dalcin   char           *cbuf;
395de78e4feSLisandro Dalcin   GmshWorkBuffer work;
396de78e4feSLisandro Dalcin   PetscErrorCode ierr;
397de78e4feSLisandro Dalcin 
398de78e4feSLisandro Dalcin   PetscFunctionBegin;
399de78e4feSLisandro Dalcin   ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr);
400de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
401de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
402de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
403de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
404de78e4feSLisandro Dalcin   ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr);
405de78e4feSLisandro Dalcin   *numVertices = (int)numTotalNodes;
406de78e4feSLisandro Dalcin   *gmsh_nodes = coordinates;
407de78e4feSLisandro Dalcin   for (v = 0, block = 0; block < numEntityBlocks; ++block) {
408de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
409de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
410de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
411de78e4feSLisandro Dalcin     if (binary) {
412de78e4feSLisandro Dalcin       int nbytes = sizeof(int) + 3*sizeof(double);
413de78e4feSLisandro Dalcin       ierr = GmshWorkBufferGet(&work, numNodes, nbytes, &cbuf);CHKERRQ(ierr);
414de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr);
415de78e4feSLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++v) {
416de78e4feSLisandro Dalcin         char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int);
417de78e4feSLisandro Dalcin         double *xyz = coordinates + v*3;
418de78e4feSLisandro Dalcin #if !defined(PETSC_WORDS_BIGENDIAN)
419de78e4feSLisandro Dalcin         ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);
420de78e4feSLisandro Dalcin         ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);
421de78e4feSLisandro Dalcin #endif
422de78e4feSLisandro Dalcin         ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr);
423de78e4feSLisandro Dalcin         ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr);
424de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
425de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
426de78e4feSLisandro Dalcin         if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift);
427de78e4feSLisandro Dalcin       }
428de78e4feSLisandro Dalcin     } else {
429de78e4feSLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++v) {
430de78e4feSLisandro Dalcin         double *xyz = coordinates + v*3;
431de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
432de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
433de78e4feSLisandro Dalcin         if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift);
434de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
435de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
436de78e4feSLisandro Dalcin       }
437de78e4feSLisandro Dalcin     }
438de78e4feSLisandro Dalcin   }
439de78e4feSLisandro Dalcin   ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr);
440de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
441de78e4feSLisandro Dalcin }
442de78e4feSLisandro Dalcin 
443de78e4feSLisandro Dalcin /*
444de78e4feSLisandro Dalcin $Elements
445de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numElements(unsigned long)
446de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long)
447de78e4feSLisandro Dalcin     tag(int) numVert[...](int)
448de78e4feSLisandro Dalcin     ...
449de78e4feSLisandro Dalcin   ...
450de78e4feSLisandro Dalcin $EndElements
451de78e4feSLisandro Dalcin */
452*730356f6SLisandro Dalcin static PetscErrorCode DMPlexCreateGmsh_ReadElements_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, GmshEntities *entities, int *numCells, GmshElement **gmsh_elems)
453de78e4feSLisandro Dalcin {
454de78e4feSLisandro Dalcin   long           c, block, numEntityBlocks, numTotalElements, elem, numElements;
455de78e4feSLisandro Dalcin   int            p, info[3], *ibuf = NULL;
456de78e4feSLisandro Dalcin   int            eid, dim, numTags, *tags, cellType, numNodes;
457*730356f6SLisandro Dalcin   GmshEntity     *entity;
458de78e4feSLisandro Dalcin   GmshElement    *elements;
459de78e4feSLisandro Dalcin   GmshWorkBuffer work;
460de78e4feSLisandro Dalcin   PetscErrorCode ierr;
461de78e4feSLisandro Dalcin 
462de78e4feSLisandro Dalcin   PetscFunctionBegin;
463de78e4feSLisandro Dalcin   ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr);
464de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
465de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
466de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
467de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);}
468de78e4feSLisandro Dalcin   ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr);
469de78e4feSLisandro Dalcin   *numCells = (int)numTotalElements;
470de78e4feSLisandro Dalcin   *gmsh_elems = elements;
471de78e4feSLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
472de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
473de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);}
474de78e4feSLisandro Dalcin     eid = info[0]; dim = info[1]; cellType = info[2];
475*730356f6SLisandro Dalcin     ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr);
476*730356f6SLisandro Dalcin     numTags = entity->numTags;
477*730356f6SLisandro Dalcin     tags = entity->tags;
478de78e4feSLisandro Dalcin     switch (cellType) {
479de78e4feSLisandro Dalcin     case 1: /* 2-node line */
480de78e4feSLisandro Dalcin       numNodes = 2;
481de78e4feSLisandro Dalcin       break;
482de78e4feSLisandro Dalcin     case 2: /* 3-node triangle */
483de78e4feSLisandro Dalcin       numNodes = 3; break;
484de78e4feSLisandro Dalcin     case 3: /* 4-node quadrangle */
485de78e4feSLisandro Dalcin       numNodes = 4;
486de78e4feSLisandro Dalcin       break;
487de78e4feSLisandro Dalcin     case 4: /* 4-node tetrahedron */
488de78e4feSLisandro Dalcin       numNodes = 4;
489de78e4feSLisandro Dalcin       break;
490de78e4feSLisandro Dalcin     case 5: /* 8-node hexahedron */
491de78e4feSLisandro Dalcin       numNodes = 8;
492de78e4feSLisandro Dalcin       break;
493de78e4feSLisandro Dalcin     case 6: /* 6-node wedge */
494de78e4feSLisandro Dalcin       numNodes = 6;
495de78e4feSLisandro Dalcin       break;
496de78e4feSLisandro Dalcin     case 15: /* 1-node vertex */
497de78e4feSLisandro Dalcin       numNodes = 1;
498de78e4feSLisandro Dalcin       break;
499de78e4feSLisandro Dalcin     default:
500de78e4feSLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
501de78e4feSLisandro Dalcin     }
502de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
503de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);}
504de78e4feSLisandro Dalcin     ierr = GmshWorkBufferGet(&work, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr);
505de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr);
506de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);}
507de78e4feSLisandro Dalcin     for (elem = 0; elem < numElements; ++elem, ++c) {
508de78e4feSLisandro Dalcin       int *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
509de78e4feSLisandro Dalcin       GmshElement *element = elements + c;
510de78e4feSLisandro Dalcin       element->dim = dim;
511de78e4feSLisandro Dalcin       element->cellType = cellType;
512de78e4feSLisandro Dalcin       element->numNodes = numNodes;
513de78e4feSLisandro Dalcin       element->numTags = numTags;
514de78e4feSLisandro Dalcin       element->id = id[0];
515de78e4feSLisandro Dalcin       for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p];
516de78e4feSLisandro Dalcin       for (p = 0; p < numTags; p++) element->tags[p] = tags[p];
517de78e4feSLisandro Dalcin     }
518de78e4feSLisandro Dalcin   }
519de78e4feSLisandro Dalcin   ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr);
520de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
521de78e4feSLisandro Dalcin }
522de78e4feSLisandro Dalcin 
523331830f3SMatthew G. Knepley /*@
5247d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
525331830f3SMatthew G. Knepley 
526331830f3SMatthew G. Knepley   Collective on comm
527331830f3SMatthew G. Knepley 
528331830f3SMatthew G. Knepley   Input Parameters:
529331830f3SMatthew G. Knepley + comm  - The MPI communicator
530331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
531331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
532331830f3SMatthew G. Knepley 
533331830f3SMatthew G. Knepley   Output Parameter:
534331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
535331830f3SMatthew G. Knepley 
536331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
5373d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
538331830f3SMatthew G. Knepley 
539331830f3SMatthew G. Knepley   Level: beginner
540331830f3SMatthew G. Knepley 
541331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
542331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
543331830f3SMatthew G. Knepley @*/
544331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
545331830f3SMatthew G. Knepley {
54611c56cb0SLisandro Dalcin   PetscViewer    parentviewer = NULL;
54711c56cb0SLisandro Dalcin   double        *coordsIn = NULL;
548*730356f6SLisandro Dalcin   GmshEntities  *entities = NULL;
5493c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
550331830f3SMatthew G. Knepley   PetscSection   coordSection;
551331830f3SMatthew G. Knepley   Vec            coordinates;
5526fbe17bfSStefano Zampini   PetscBT        periodicV = NULL, periodicC = NULL;
55384572febSToby Isaac   PetscScalar   *coords;
554dea2a3c7SStefano Zampini   PetscInt       dim = 0, embedDim = -1, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE;
5551d64f2ccSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1;
55611c56cb0SLisandro Dalcin   PetscMPIInt    rank;
557331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
558dea2a3c7SStefano Zampini   PetscBool      binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE;
559dea2a3c7SStefano Zampini   PetscBool      enable_hybrid = PETSC_FALSE;
560331830f3SMatthew G. Knepley   PetscErrorCode ierr;
561331830f3SMatthew G. Knepley 
562331830f3SMatthew G. Knepley   PetscFunctionBegin;
563331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
564dea2a3c7SStefano Zampini   ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr);
565dea2a3c7SStefano Zampini   ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr);
566dea2a3c7SStefano Zampini   ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr);
567dea2a3c7SStefano Zampini   ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr);
568dea2a3c7SStefano Zampini   ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr);
56911c56cb0SLisandro Dalcin   if (zerobase) shift = 0;
57011c56cb0SLisandro Dalcin 
571331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
572331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
5733b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
57411c56cb0SLisandro Dalcin 
57511c56cb0SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
57611c56cb0SLisandro Dalcin 
57711c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
5783b46f529SLisandro Dalcin   if (binary) {
57911c56cb0SLisandro Dalcin     parentviewer = viewer;
58011c56cb0SLisandro Dalcin     ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
58111c56cb0SLisandro Dalcin   }
58211c56cb0SLisandro Dalcin 
5833b46f529SLisandro Dalcin   if (!rank) {
584dea2a3c7SStefano Zampini     PetscBool match, hybrid;
585de78e4feSLisandro Dalcin     int       fileType, fileFormat, dataSize;
586f6ac7a6aSMichael Lange     float     version;
587331830f3SMatthew G. Knepley 
588331830f3SMatthew G. Knepley     /* Read header */
589060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
59004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
591331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
592060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
593f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
594f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
595f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
596de78e4feSLisandro Dalcin     if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported");
597de78e4feSLisandro Dalcin     if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0");
598331830f3SMatthew 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);
599de78e4feSLisandro Dalcin     fileFormat = (int)version;
60004d1ad83SMichael Lange     if (binary) {
601b9eae255SMichael Lange       int checkInt;
602060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
60304d1ad83SMichael Lange       if (checkInt != 1) {
604b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
60511c56cb0SLisandro Dalcin         if (checkInt == 1) byteSwap = PETSC_TRUE;
60604d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
60704d1ad83SMichael Lange       }
6080877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
609060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
61004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
611331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
61211c56cb0SLisandro Dalcin 
613dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
614060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
615dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
616dc0ede3bSMatthew G. Knepley     if (match) {
617dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
618dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
619dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
62011c56cb0SLisandro Dalcin       for (r = 0; r < numRegions; ++r) {
62111c56cb0SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
62211c56cb0SLisandro Dalcin       }
623dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
624dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
625dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
626dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
627dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
628dc0ede3bSMatthew G. Knepley     }
62911c56cb0SLisandro Dalcin 
630de78e4feSLisandro Dalcin     /* Read entities */
631de78e4feSLisandro Dalcin     if (fileFormat == 4) {
632de78e4feSLisandro Dalcin       ierr = PetscStrncmp(line, "$Entities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
633de78e4feSLisandro Dalcin       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
634*730356f6SLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadEntities_v4(viewer, binary, byteSwap, &entities);CHKERRQ(ierr);
635de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
636de78e4feSLisandro Dalcin       ierr = PetscStrncmp(line, "$EndEntities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
637de78e4feSLisandro Dalcin       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
638de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
639de78e4feSLisandro Dalcin     }
640de78e4feSLisandro Dalcin 
641dc0ede3bSMatthew G. Knepley     /* Read vertices */
64204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
643331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
644de78e4feSLisandro Dalcin     if (fileFormat == 2) {
645de78e4feSLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadNodes_v2(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr);
646de78e4feSLisandro Dalcin     } else {
647de78e4feSLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadNodes_v4(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr);
648de78e4feSLisandro Dalcin     }
649060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
65004d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
651331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
65211c56cb0SLisandro Dalcin 
653331830f3SMatthew G. Knepley     /* Read cells */
654060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
65504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
656331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
657de78e4feSLisandro Dalcin     if (fileFormat == 2) {
658de78e4feSLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadElements_v2(viewer, binary, byteSwap, shift, &numCells, &gmsh_elem);CHKERRQ(ierr);
659de78e4feSLisandro Dalcin     } else {
660de78e4feSLisandro Dalcin       ierr = DMPlexCreateGmsh_ReadElements_v4(viewer, binary, byteSwap, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr);
661de78e4feSLisandro Dalcin     }
662de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
663de78e4feSLisandro Dalcin     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
664de78e4feSLisandro Dalcin     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
665de78e4feSLisandro Dalcin 
666*730356f6SLisandro Dalcin     ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr);
667de78e4feSLisandro Dalcin 
668dea2a3c7SStefano Zampini     hybrid = PETSC_FALSE;
669a4bb7517SMichael Lange     for (trueNumCells = 0, c = 0; c < numCells; ++c) {
670dea2a3c7SStefano Zampini       int on = -1;
671a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
672dea2a3c7SStefano Zampini       if (gmsh_elem[c].dim == dim) {hybrid = (trueNumCells ? (on != gmsh_elem[c].cellType ? on = gmsh_elem[c].cellType,PETSC_TRUE : hybrid) : (on = gmsh_elem[c].cellType, PETSC_FALSE) ); trueNumCells++;}
673dea2a3c7SStefano Zampini       /* wedges always indicate an hybrid mesh in PLEX */
674dea2a3c7SStefano Zampini       if (on == 6 || on == 13) hybrid = PETSC_TRUE;
675db397164SMichael Lange     }
67611c56cb0SLisandro Dalcin 
677dea2a3c7SStefano Zampini     /* Renumber cells for hybrid grids */
678dea2a3c7SStefano Zampini     if (hybrid && enable_hybrid) {
679dea2a3c7SStefano Zampini       PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL;
680dea2a3c7SStefano Zampini       PetscInt cell, tn, *tp;
681dea2a3c7SStefano Zampini       int      n1 = 0,n2 = 0;
682dea2a3c7SStefano Zampini 
683dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr);
684dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr);
685dea2a3c7SStefano Zampini       for (cell = 0, c = 0; c < numCells; ++c) {
686dea2a3c7SStefano Zampini         if (gmsh_elem[c].dim == dim) {
687dea2a3c7SStefano Zampini           if (!n1) n1 = gmsh_elem[c].cellType;
688dea2a3c7SStefano Zampini           else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType;
689dea2a3c7SStefano Zampini 
690dea2a3c7SStefano Zampini           if      (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell;
691dea2a3c7SStefano Zampini           else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell;
692dea2a3c7SStefano Zampini           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types");
693dea2a3c7SStefano Zampini           cell++;
694dea2a3c7SStefano Zampini         }
695dea2a3c7SStefano Zampini       }
696dea2a3c7SStefano Zampini 
697dea2a3c7SStefano Zampini       switch (n1) {
698dea2a3c7SStefano Zampini       case 2: /* triangles */
699dea2a3c7SStefano Zampini       case 9:
700dea2a3c7SStefano Zampini         switch (n2) {
701dea2a3c7SStefano Zampini         case 0: /* single type mesh */
702dea2a3c7SStefano Zampini         case 3: /* quads */
703dea2a3c7SStefano Zampini         case 10:
704dea2a3c7SStefano Zampini           break;
705dea2a3c7SStefano Zampini         default:
706dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
707dea2a3c7SStefano Zampini         }
708dea2a3c7SStefano Zampini         break;
709dea2a3c7SStefano Zampini       case 3:
710dea2a3c7SStefano Zampini       case 10:
711dea2a3c7SStefano Zampini         switch (n2) {
712dea2a3c7SStefano Zampini         case 0: /* single type mesh */
713dea2a3c7SStefano Zampini         case 2: /* swap since we list simplices first */
714dea2a3c7SStefano Zampini         case 9:
715dea2a3c7SStefano Zampini           tn  = hc1;
716dea2a3c7SStefano Zampini           hc1 = hc2;
717dea2a3c7SStefano Zampini           hc2 = tn;
718dea2a3c7SStefano Zampini 
719dea2a3c7SStefano Zampini           tp           = hybridCells1;
720dea2a3c7SStefano Zampini           hybridCells1 = hybridCells2;
721dea2a3c7SStefano Zampini           hybridCells2 = tp;
722dea2a3c7SStefano Zampini           break;
723dea2a3c7SStefano Zampini         default:
724dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
725dea2a3c7SStefano Zampini         }
726dea2a3c7SStefano Zampini         break;
727dea2a3c7SStefano Zampini       case 4: /* tetrahedra */
728dea2a3c7SStefano Zampini       case 11:
729dea2a3c7SStefano Zampini         switch (n2) {
730dea2a3c7SStefano Zampini         case 0: /* single type mesh */
731dea2a3c7SStefano Zampini         case 6: /* wedges */
732dea2a3c7SStefano Zampini         case 13:
733dea2a3c7SStefano Zampini           break;
734dea2a3c7SStefano Zampini         default:
735dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
736dea2a3c7SStefano Zampini         }
737dea2a3c7SStefano Zampini         break;
738dea2a3c7SStefano Zampini       case 6:
739dea2a3c7SStefano Zampini       case 13:
740dea2a3c7SStefano Zampini         switch (n2) {
741dea2a3c7SStefano Zampini         case 0: /* single type mesh */
742dea2a3c7SStefano Zampini         case 4: /* swap since we list simplices first */
743dea2a3c7SStefano Zampini         case 11:
744dea2a3c7SStefano Zampini           tn  = hc1;
745dea2a3c7SStefano Zampini           hc1 = hc2;
746dea2a3c7SStefano Zampini           hc2 = tn;
747dea2a3c7SStefano Zampini 
748dea2a3c7SStefano Zampini           tp           = hybridCells1;
749dea2a3c7SStefano Zampini           hybridCells1 = hybridCells2;
750dea2a3c7SStefano Zampini           hybridCells2 = tp;
751dea2a3c7SStefano Zampini           break;
752dea2a3c7SStefano Zampini         default:
753dea2a3c7SStefano Zampini           SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
754dea2a3c7SStefano Zampini         }
755dea2a3c7SStefano Zampini         break;
756dea2a3c7SStefano Zampini       default:
757dea2a3c7SStefano Zampini         SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2);
758dea2a3c7SStefano Zampini       }
759dea2a3c7SStefano Zampini       cMax = hc1;
760dea2a3c7SStefano Zampini       ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr);
761dea2a3c7SStefano Zampini       for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell;
762dea2a3c7SStefano Zampini       for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1;
763dea2a3c7SStefano Zampini       ierr = PetscFree(hybridCells1);CHKERRQ(ierr);
764dea2a3c7SStefano Zampini       ierr = PetscFree(hybridCells2);CHKERRQ(ierr);
765dea2a3c7SStefano Zampini     }
766dea2a3c7SStefano Zampini 
76711c56cb0SLisandro Dalcin     /* OPTIONAL Read periodic section */
768d08df55aSStefano Zampini     if (periodic) {
769f45c9589SStefano Zampini       PetscInt pVert, *periodicMapT, *aux;
770d08df55aSStefano Zampini       int      numPeriodic;
771d08df55aSStefano Zampini 
772f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
7736fbe17bfSStefano Zampini       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
774f45c9589SStefano Zampini       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
775de78e4feSLisandro Dalcin 
776de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
777de78e4feSLisandro Dalcin       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
778de78e4feSLisandro Dalcin       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
779de78e4feSLisandro Dalcin       if (fileFormat == 2 || !binary) {
780d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
781d08df55aSStefano Zampini         snum = sscanf(line, "%d", &numPeriodic);
782d08df55aSStefano Zampini         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
783de78e4feSLisandro Dalcin       } else {
784de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
785de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);}
786de78e4feSLisandro Dalcin       }
787d08df55aSStefano Zampini       for (i = 0; i < numPeriodic; i++) {
788de78e4feSLisandro Dalcin         int    ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode;
789de78e4feSLisandro Dalcin         long   j, nNodes;
790de78e4feSLisandro Dalcin         double affine[16];
791d08df55aSStefano Zampini 
792de78e4feSLisandro Dalcin         if (fileFormat == 2 || !binary) {
793d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
794de78e4feSLisandro Dalcin           snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag);
795d08df55aSStefano Zampini           if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
796de78e4feSLisandro Dalcin         } else {
797de78e4feSLisandro Dalcin           ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
798de78e4feSLisandro Dalcin           if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
799de78e4feSLisandro Dalcin           slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2];
800de78e4feSLisandro Dalcin         }
801de78e4feSLisandro Dalcin         (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */
802de78e4feSLisandro Dalcin 
803de78e4feSLisandro Dalcin         if (fileFormat == 2 || !binary) {
804da83f57bSLisandro Dalcin           ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
805de78e4feSLisandro Dalcin           snum = sscanf(line, "%ld", &nNodes);
806da83f57bSLisandro Dalcin           if (snum != 1) { /* discard tranformation and try again */
80772ffbcc9SStefano Zampini             ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
808d08df55aSStefano Zampini             ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
809de78e4feSLisandro Dalcin             snum = sscanf(line, "%ld", &nNodes);
810d08df55aSStefano Zampini             if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
811da83f57bSLisandro Dalcin           }
812de78e4feSLisandro Dalcin         } else {
813de78e4feSLisandro Dalcin           ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
814de78e4feSLisandro Dalcin           if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
815de78e4feSLisandro Dalcin           if (nNodes == -1) { /* discard tranformation and try again */
816de78e4feSLisandro Dalcin             ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
817de78e4feSLisandro Dalcin             ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
818de78e4feSLisandro Dalcin             if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
819de78e4feSLisandro Dalcin           }
820de78e4feSLisandro Dalcin         }
821de78e4feSLisandro Dalcin 
822d08df55aSStefano Zampini         for (j = 0; j < nNodes; j++) {
823de78e4feSLisandro Dalcin           if (fileFormat == 2 || !binary) {
824d08df55aSStefano Zampini             ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
82511c56cb0SLisandro Dalcin             snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
826d08df55aSStefano Zampini             if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
827de78e4feSLisandro Dalcin           } else {
828de78e4feSLisandro Dalcin             ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr);
829de78e4feSLisandro Dalcin             if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);}
830de78e4feSLisandro Dalcin             slaveNode = ibuf[0]; masterNode = ibuf[1];
831de78e4feSLisandro Dalcin           }
832917266a4SLisandro Dalcin           periodicMapT[slaveNode - shift] = masterNode - shift;
833917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr);
834917266a4SLisandro Dalcin           ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr);
835d08df55aSStefano Zampini         }
836d08df55aSStefano Zampini       }
837d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
838d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
839d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
840de78e4feSLisandro Dalcin 
841d08df55aSStefano Zampini       /* we may have slaves of slaves */
842d08df55aSStefano Zampini       for (i = 0; i < numVertices; i++) {
843f45c9589SStefano Zampini         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
844f45c9589SStefano Zampini           periodicMapT[i] = periodicMapT[periodicMapT[i]];
845d08df55aSStefano Zampini         }
846d08df55aSStefano Zampini       }
847f45c9589SStefano Zampini       /* periodicMap : from old to new numbering (periodic vertices excluded)
848f45c9589SStefano Zampini          periodicMapI: from new to old numbering */
849f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
850f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
851f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
852f45c9589SStefano Zampini       for (i = 0, pVert = 0; i < numVertices; i++) {
853f45c9589SStefano Zampini         if (periodicMapT[i] != i) {
854f45c9589SStefano Zampini           pVert++;
855f45c9589SStefano Zampini         } else {
856f45c9589SStefano Zampini           aux[i] = i - pVert;
857f45c9589SStefano Zampini           periodicMapI[i - pVert] = i;
858f45c9589SStefano Zampini         }
859f45c9589SStefano Zampini       }
860f45c9589SStefano Zampini       for (i = 0 ; i < numVertices; i++) {
861f45c9589SStefano Zampini         periodicMap[i] = aux[periodicMapT[i]];
862f45c9589SStefano Zampini       }
863f45c9589SStefano Zampini       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
864f45c9589SStefano Zampini       ierr = PetscFree(aux);CHKERRQ(ierr);
865f45c9589SStefano Zampini       /* remove periodic vertices */
866f45c9589SStefano Zampini       numVertices = numVertices - pVert;
867d08df55aSStefano Zampini     }
868db397164SMichael Lange   }
86911c56cb0SLisandro Dalcin 
87011c56cb0SLisandro Dalcin   if (parentviewer) {
87111c56cb0SLisandro Dalcin     ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
87211c56cb0SLisandro Dalcin   }
87311c56cb0SLisandro Dalcin 
874a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
875db397164SMichael Lange   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
876a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
877a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
878dea2a3c7SStefano Zampini       ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
879a4bb7517SMichael Lange       cell++;
880db397164SMichael Lange     }
881331830f3SMatthew G. Knepley   }
882331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
883a4bb7517SMichael Lange   /* Add cell-vertex connections */
884a4bb7517SMichael Lange   for (cell = 0, c = 0; c < numCells; ++c) {
885a4bb7517SMichael Lange     if (gmsh_elem[c].dim == dim) {
88611c56cb0SLisandro Dalcin       PetscInt pcone[8], corner;
887a4bb7517SMichael Lange       for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
888dd169d64SMatthew G. Knepley         const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
889917266a4SLisandro Dalcin         pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells;
890db397164SMichael Lange       }
89197ed6be6Ssarens       if (dim == 3) {
89297ed6be6Ssarens         /* Tetrahedra are inverted */
893dea2a3c7SStefano Zampini         if (gmsh_elem[c].cellType == 4) {
89497ed6be6Ssarens           PetscInt tmp = pcone[0];
89597ed6be6Ssarens           pcone[0] = pcone[1];
89697ed6be6Ssarens           pcone[1] = tmp;
89797ed6be6Ssarens         }
89897ed6be6Ssarens         /* Hexahedra are inverted */
899dea2a3c7SStefano Zampini         if (gmsh_elem[c].cellType == 5) {
90097ed6be6Ssarens           PetscInt tmp = pcone[1];
90197ed6be6Ssarens           pcone[1] = pcone[3];
90297ed6be6Ssarens           pcone[3] = tmp;
90397ed6be6Ssarens         }
904dea2a3c7SStefano Zampini         /* Prisms are inverted */
905dea2a3c7SStefano Zampini         if (gmsh_elem[c].cellType == 6) {
906dea2a3c7SStefano Zampini           PetscInt tmp;
907dea2a3c7SStefano Zampini 
908dea2a3c7SStefano Zampini           tmp      = pcone[1];
909dea2a3c7SStefano Zampini           pcone[1] = pcone[2];
910dea2a3c7SStefano Zampini           pcone[2] = tmp;
911dea2a3c7SStefano Zampini           tmp      = pcone[4];
912dea2a3c7SStefano Zampini           pcone[4] = pcone[5];
913dea2a3c7SStefano Zampini           pcone[5] = tmp;
91497ed6be6Ssarens         }
915dea2a3c7SStefano Zampini       } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */
916dea2a3c7SStefano Zampini         /* quads are input to PLEX as prisms */
917dea2a3c7SStefano Zampini         if (gmsh_elem[c].cellType == 3) {
918dea2a3c7SStefano Zampini           PetscInt tmp = pcone[2];
919dea2a3c7SStefano Zampini           pcone[2] = pcone[3];
920dea2a3c7SStefano Zampini           pcone[3] = tmp;
921dea2a3c7SStefano Zampini         }
922dea2a3c7SStefano Zampini       }
923dea2a3c7SStefano Zampini       ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr);
924a4bb7517SMichael Lange       cell++;
925331830f3SMatthew G. Knepley     }
926a4bb7517SMichael Lange   }
9276228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
928c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
929dea2a3c7SStefano Zampini   ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
930331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
931331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
932331830f3SMatthew G. Knepley   if (interpolate) {
9335fd9971aSMatthew G. Knepley     DM idm;
934331830f3SMatthew G. Knepley 
935331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
936331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
937331830f3SMatthew G. Knepley     *dm  = idm;
938331830f3SMatthew G. Knepley   }
939917266a4SLisandro Dalcin 
940f6c8a31fSLisandro Dalcin   if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
941917266a4SLisandro Dalcin   if (!rank && usemarker) {
942d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
943d3f73514SStefano Zampini 
944d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
945d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
946d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
947d3f73514SStefano Zampini       PetscInt suppSize;
948d3f73514SStefano Zampini 
949d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
950d3f73514SStefano Zampini       if (suppSize == 1) {
951d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
952d3f73514SStefano Zampini 
953d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
954d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
955d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
956d3f73514SStefano Zampini         }
957d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
958d3f73514SStefano Zampini       }
959d3f73514SStefano Zampini     }
960d3f73514SStefano Zampini   }
96116ad7e67SMichael Lange 
96216ad7e67SMichael Lange   if (!rank) {
96311c56cb0SLisandro Dalcin     PetscInt vStart, vEnd;
964d1a54cefSMatthew G. Knepley 
96516ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
96611c56cb0SLisandro Dalcin     for (cell = 0, c = 0; c < numCells; ++c) {
96711c56cb0SLisandro Dalcin 
96811c56cb0SLisandro Dalcin       /* Create face sets */
9695964b3dcSLisandro Dalcin       if (interpolate && gmsh_elem[c].dim == dim-1) {
97016ad7e67SMichael Lange         const PetscInt *join;
97111c56cb0SLisandro Dalcin         PetscInt        joinSize, pcone[8], corner;
97211c56cb0SLisandro Dalcin         /* Find the relevant facet with vertex joins */
973a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
974dd169d64SMatthew G. Knepley           const PetscInt cc = gmsh_elem[c].nodes[corner] - shift;
975917266a4SLisandro Dalcin           pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart;
97616ad7e67SMichael Lange         }
97711c56cb0SLisandro Dalcin         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr);
978f10f1c67SMatthew G. Knepley         if (joinSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for GMsh element %d (Plex cell %D)", gmsh_elem[c].id, c);
979c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
980a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
98116ad7e67SMichael Lange       }
9826e1dd89aSLawrence Mitchell 
9836e1dd89aSLawrence Mitchell       /* Create cell sets */
9846e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
9856e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
986dea2a3c7SStefano Zampini           ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
9876e1dd89aSLawrence Mitchell         }
988917266a4SLisandro Dalcin         cell++;
9896e1dd89aSLawrence Mitchell       }
9900c070f12SSander Arens 
9910c070f12SSander Arens       /* Create vertex sets */
9920c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
9930c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
994917266a4SLisandro Dalcin           const PetscInt cc = gmsh_elem[c].nodes[0] - shift;
99511c56cb0SLisandro Dalcin           const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart;
996d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
9970c070f12SSander Arens         }
9980c070f12SSander Arens       }
9990c070f12SSander Arens     }
100016ad7e67SMichael Lange   }
100116ad7e67SMichael Lange 
100211c56cb0SLisandro Dalcin   /* Create coordinates */
1003dea2a3c7SStefano Zampini   if (embedDim < 0) embedDim = dim;
1004dea2a3c7SStefano Zampini   ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr);
1005979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1006331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
10071d64f2ccSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr);
1008f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
1009f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
1010f45c9589SStefano Zampini   } else {
101175b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
1012f45c9589SStefano Zampini   }
101375b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
10141d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr);
10151d64f2ccSMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr);
1016331830f3SMatthew G. Knepley   }
101711c56cb0SLisandro Dalcin   if (periodicMap) {
1018437bdd5bSStefano Zampini     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
1019f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
1020f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
1021437bdd5bSStefano Zampini         PetscInt  corner;
102211c56cb0SLisandro Dalcin         PetscBool pc = PETSC_FALSE;
1023437bdd5bSStefano Zampini         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1024917266a4SLisandro Dalcin           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift));
1025437bdd5bSStefano Zampini         }
1026437bdd5bSStefano Zampini         if (pc) {
102711c56cb0SLisandro Dalcin           PetscInt dof = gmsh_elem[c].numNodes*embedDim;
1028dea2a3c7SStefano Zampini           PetscInt ucell = hybridMap ? hybridMap[cell] : cell;
1029dea2a3c7SStefano Zampini           ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr);
1030dea2a3c7SStefano Zampini           ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr);
1031dea2a3c7SStefano Zampini           ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr);
10326fbe17bfSStefano Zampini         }
1033f45c9589SStefano Zampini         cell++;
1034f45c9589SStefano Zampini       }
1035f45c9589SStefano Zampini     }
1036f45c9589SStefano Zampini   }
1037331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1038331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
10398b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1040331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1041331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
10421d64f2ccSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr);
1043331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
1044331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
1045f45c9589SStefano Zampini   if (periodicMap) {
1046f45c9589SStefano Zampini     PetscInt off;
1047f45c9589SStefano Zampini 
1048f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
1049f45c9589SStefano Zampini       PetscInt pcone[8], corner;
1050f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
1051dea2a3c7SStefano Zampini         PetscInt ucell = hybridMap ? hybridMap[cell] : cell;
1052dea2a3c7SStefano Zampini         if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) {
1053f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
1054917266a4SLisandro Dalcin             pcone[corner] = gmsh_elem[c].nodes[corner] - shift;
1055f45c9589SStefano Zampini           }
1056f45c9589SStefano Zampini           if (dim == 3) {
1057f45c9589SStefano Zampini             /* Tetrahedra are inverted */
1058dea2a3c7SStefano Zampini             if (gmsh_elem[c].cellType == 4) {
1059f45c9589SStefano Zampini               PetscInt tmp = pcone[0];
1060f45c9589SStefano Zampini               pcone[0] = pcone[1];
1061f45c9589SStefano Zampini               pcone[1] = tmp;
1062f45c9589SStefano Zampini             }
1063f45c9589SStefano Zampini             /* Hexahedra are inverted */
1064dea2a3c7SStefano Zampini             if (gmsh_elem[c].cellType == 5) {
1065f45c9589SStefano Zampini               PetscInt tmp = pcone[1];
1066f45c9589SStefano Zampini               pcone[1] = pcone[3];
1067f45c9589SStefano Zampini               pcone[3] = tmp;
1068f45c9589SStefano Zampini             }
1069dea2a3c7SStefano Zampini             /* Prisms are inverted */
1070dea2a3c7SStefano Zampini             if (gmsh_elem[c].cellType == 6) {
1071dea2a3c7SStefano Zampini               PetscInt tmp;
1072dea2a3c7SStefano Zampini 
1073dea2a3c7SStefano Zampini               tmp      = pcone[1];
1074dea2a3c7SStefano Zampini               pcone[1] = pcone[2];
1075dea2a3c7SStefano Zampini               pcone[2] = tmp;
1076dea2a3c7SStefano Zampini               tmp      = pcone[4];
1077dea2a3c7SStefano Zampini               pcone[4] = pcone[5];
1078dea2a3c7SStefano Zampini               pcone[5] = tmp;
1079f45c9589SStefano Zampini             }
1080dea2a3c7SStefano Zampini           } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */
1081dea2a3c7SStefano Zampini             /* quads are input to PLEX as prisms */
1082dea2a3c7SStefano Zampini             if (gmsh_elem[c].cellType == 3) {
1083dea2a3c7SStefano Zampini               PetscInt tmp = pcone[2];
1084dea2a3c7SStefano Zampini               pcone[2] = pcone[3];
1085dea2a3c7SStefano Zampini               pcone[3] = tmp;
1086dea2a3c7SStefano Zampini             }
1087dea2a3c7SStefano Zampini           }
1088dea2a3c7SStefano Zampini           ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr);
1089f45c9589SStefano Zampini           for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
109011c56cb0SLisandro Dalcin             v = pcone[corner];
1091dd169d64SMatthew G. Knepley             for (d = 0; d < embedDim; ++d) {
109211c56cb0SLisandro Dalcin               coords[off++] = (PetscReal) coordsIn[v*3+d];
1093f45c9589SStefano Zampini             }
1094f45c9589SStefano Zampini           }
10956fbe17bfSStefano Zampini         }
1096f45c9589SStefano Zampini         cell++;
1097f45c9589SStefano Zampini       }
1098f45c9589SStefano Zampini     }
1099f45c9589SStefano Zampini     for (v = 0; v < numVertices; ++v) {
1100f45c9589SStefano Zampini       ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
1101dd169d64SMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
110211c56cb0SLisandro Dalcin         coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d];
1103f45c9589SStefano Zampini       }
1104f45c9589SStefano Zampini     }
1105f45c9589SStefano Zampini   } else {
1106331830f3SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
11071d64f2ccSMatthew G. Knepley       for (d = 0; d < embedDim; ++d) {
110811c56cb0SLisandro Dalcin         coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d];
1109331830f3SMatthew G. Knepley       }
1110331830f3SMatthew G. Knepley     }
1111331830f3SMatthew G. Knepley   }
1112331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1113331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
111411c56cb0SLisandro Dalcin   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
111511c56cb0SLisandro Dalcin 
111611c56cb0SLisandro Dalcin   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
111711c56cb0SLisandro Dalcin   ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
1118dea2a3c7SStefano Zampini   ierr = PetscFree(hybridMap);CHKERRQ(ierr);
1119d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
1120fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
11216fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
11226fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
112311c56cb0SLisandro Dalcin   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
112411c56cb0SLisandro Dalcin 
11253b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
1126331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
1127331830f3SMatthew G. Knepley }
1128