xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 3548e9b60dcb5707af254b16ca03bdb62b941e5c)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3331830f3SMatthew G. Knepley 
47d282ae0SMichael Lange /*@C
57d282ae0SMichael Lange   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
67d282ae0SMichael Lange 
77d282ae0SMichael Lange + comm        - The MPI communicator
87d282ae0SMichael Lange . filename    - Name of the Gmsh file
97d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh
107d282ae0SMichael Lange 
117d282ae0SMichael Lange   Output Parameter:
127d282ae0SMichael Lange . dm  - The DM object representing the mesh
137d282ae0SMichael Lange 
147d282ae0SMichael Lange   Level: beginner
157d282ae0SMichael Lange 
167d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
177d282ae0SMichael Lange @*/
187d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
197d282ae0SMichael Lange {
207d282ae0SMichael Lange   PetscViewer     viewer, vheader;
21abc86ac4SMichael Lange   PetscMPIInt     rank;
227d282ae0SMichael Lange   PetscViewerType vtype;
237d282ae0SMichael Lange   char            line[PETSC_MAX_PATH_LEN];
247d282ae0SMichael Lange   int             snum;
257d282ae0SMichael Lange   PetscBool       match;
26c1c22fd2SMatthew G. Knepley   int             fT;
27c1c22fd2SMatthew G. Knepley   PetscInt        fileType;
28f6ac7a6aSMichael Lange   float           version;
297d282ae0SMichael Lange   PetscErrorCode  ierr;
307d282ae0SMichael Lange 
317d282ae0SMichael Lange   PetscFunctionBegin;
32abc86ac4SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
337d282ae0SMichael Lange   /* Determine Gmsh file type (ASCII or binary) from file header */
347d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr);
357d282ae0SMichael Lange   ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr);
367d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr);
377d282ae0SMichael Lange   ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr);
38abc86ac4SMichael Lange   if (!rank) {
397d282ae0SMichael Lange     /* Read only the first two lines of the Gmsh file */
40060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
417d282ae0SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
427d282ae0SMichael Lange     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
43060da220SMatthew G. Knepley     ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
44c1c22fd2SMatthew G. Knepley     snum = sscanf(line, "%f %d", &version, &fT);
45c1c22fd2SMatthew G. Knepley     fileType = (PetscInt) fT;
46f6ac7a6aSMichael Lange     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
47f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
48abc86ac4SMichael Lange   }
49abc86ac4SMichael Lange   ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
507d282ae0SMichael Lange   /* Create appropriate viewer and build plex */
517d282ae0SMichael Lange   if (fileType == 0) vtype = PETSCVIEWERASCII;
527d282ae0SMichael Lange   else vtype = PETSCVIEWERBINARY;
537d282ae0SMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
547d282ae0SMichael Lange   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
557d282ae0SMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
567d282ae0SMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
577d282ae0SMichael Lange   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
587d282ae0SMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
597d282ae0SMichael Lange   ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr);
607d282ae0SMichael Lange   PetscFunctionReturn(0);
617d282ae0SMichael Lange }
627d282ae0SMichael Lange 
63df032b82SMatthew G. Knepley static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems)
64df032b82SMatthew G. Knepley {
65df032b82SMatthew G. Knepley   PetscInt       c, p;
66df032b82SMatthew G. Knepley   GmshElement   *elements;
67bf6ba3a3SMatthew G. Knepley   int            i, cellType, dim, numNodes, numNodesIgnore, numElem, numTags;
68bf6ba3a3SMatthew G. Knepley   PetscInt       pibuf[64];
69df032b82SMatthew G. Knepley   int            ibuf[16];
70df032b82SMatthew G. Knepley   PetscErrorCode ierr;
71df032b82SMatthew G. Knepley 
72df032b82SMatthew G. Knepley   PetscFunctionBegin;
73df032b82SMatthew G. Knepley   ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr);
74df032b82SMatthew G. Knepley   for (c = 0; c < numCells;) {
75df032b82SMatthew G. Knepley     ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
76df032b82SMatthew G. Knepley     if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);
77df032b82SMatthew G. Knepley     if (binary) {
78df032b82SMatthew G. Knepley       cellType = ibuf[0];
79df032b82SMatthew G. Knepley       numElem = ibuf[1];
80df032b82SMatthew G. Knepley       numTags = ibuf[2];
81df032b82SMatthew G. Knepley     } else {
82df032b82SMatthew G. Knepley       elements[c].id = ibuf[0];
83df032b82SMatthew G. Knepley       cellType = ibuf[1];
84df032b82SMatthew G. Knepley       numTags = ibuf[2];
85df032b82SMatthew G. Knepley       numElem = 1;
86df032b82SMatthew G. Knepley     }
87bf6ba3a3SMatthew G. Knepley     /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */
88bf6ba3a3SMatthew G. Knepley     numNodesIgnore = 0;
89df032b82SMatthew G. Knepley     switch (cellType) {
90df032b82SMatthew G. Knepley     case 1: /* 2-node line */
91df032b82SMatthew G. Knepley       dim = 1;
92df032b82SMatthew G. Knepley       numNodes = 2;
93df032b82SMatthew G. Knepley       break;
94df032b82SMatthew G. Knepley     case 2: /* 3-node triangle */
95df032b82SMatthew G. Knepley       dim = 2;
96df032b82SMatthew G. Knepley       numNodes = 3;
97df032b82SMatthew G. Knepley       break;
98df032b82SMatthew G. Knepley     case 3: /* 4-node quadrangle */
99df032b82SMatthew G. Knepley       dim = 2;
100df032b82SMatthew G. Knepley       numNodes = 4;
101df032b82SMatthew G. Knepley       break;
102df032b82SMatthew G. Knepley     case 4: /* 4-node tetrahedron */
103df032b82SMatthew G. Knepley       dim  = 3;
104df032b82SMatthew G. Knepley       numNodes = 4;
105df032b82SMatthew G. Knepley       break;
106df032b82SMatthew G. Knepley     case 5: /* 8-node hexahedron */
107df032b82SMatthew G. Knepley       dim = 3;
108df032b82SMatthew G. Knepley       numNodes = 8;
109df032b82SMatthew G. Knepley       break;
110bf6ba3a3SMatthew G. Knepley     case 8: /* 3-node 2nd order line */
111bf6ba3a3SMatthew G. Knepley       dim = 1;
112bf6ba3a3SMatthew G. Knepley       numNodes = 2;
113bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 1;
114bf6ba3a3SMatthew G. Knepley       break;
115bf6ba3a3SMatthew G. Knepley     case 9: /* 6-node 2nd order triangle */
116bf6ba3a3SMatthew G. Knepley       dim = 2;
117bf6ba3a3SMatthew G. Knepley       numNodes = 3;
118bf6ba3a3SMatthew G. Knepley       numNodesIgnore = 3;
119bf6ba3a3SMatthew G. Knepley       break;
120df032b82SMatthew G. Knepley     case 15: /* 1-node vertex */
121df032b82SMatthew G. Knepley       dim = 0;
122df032b82SMatthew G. Knepley       numNodes = 1;
123df032b82SMatthew G. Knepley       break;
124bf6ba3a3SMatthew G. Knepley     case 6: /* 6-node prism */
125bf6ba3a3SMatthew G. Knepley     case 7: /* 5-node pyramid */
126bf6ba3a3SMatthew G. Knepley     case 10: /* 9-node 2nd order quadrangle */
127bf6ba3a3SMatthew G. Knepley     case 11: /* 10-node 2nd order tetrahedron */
128bf6ba3a3SMatthew G. Knepley     case 12: /* 27-node 2nd order hexhedron */
129bf6ba3a3SMatthew G. Knepley     case 13: /* 19-node 2nd order prism */
130bf6ba3a3SMatthew G. Knepley     case 14: /* 14-node 2nd order pyramid */
131df032b82SMatthew G. Knepley     default:
132df032b82SMatthew G. Knepley       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType);
133df032b82SMatthew G. Knepley     }
134df032b82SMatthew G. Knepley     if (binary) {
135bf6ba3a3SMatthew G. Knepley       const PetscInt nint = numNodes + numTags + 1 + numNodesIgnore;
136df032b82SMatthew G. Knepley       for (i = 0; i < numElem; ++i, ++c) {
137df032b82SMatthew G. Knepley         /* Loop over inner binary element block */
138df032b82SMatthew G. Knepley         elements[c].dim = dim;
139df032b82SMatthew G. Knepley         elements[c].numNodes = numNodes;
140df032b82SMatthew G. Knepley         elements[c].numTags = numTags;
141df032b82SMatthew G. Knepley 
142df032b82SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
143df032b82SMatthew G. Knepley         if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);
144df032b82SMatthew G. Knepley         elements[c].id = ibuf[0];
145df032b82SMatthew G. Knepley         for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p];
146df032b82SMatthew G. Knepley         for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p];
147df032b82SMatthew G. Knepley       }
148df032b82SMatthew G. Knepley     } else {
149df032b82SMatthew G. Knepley       elements[c].dim = dim;
150df032b82SMatthew G. Knepley       elements[c].numNodes = numNodes;
151df032b82SMatthew G. Knepley       elements[c].numTags = numTags;
152df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr);
153df032b82SMatthew G. Knepley       ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr);
154bf6ba3a3SMatthew G. Knepley       ierr = PetscViewerRead(viewer, pibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr);
155df032b82SMatthew G. Knepley       c++;
156df032b82SMatthew G. Knepley     }
157df032b82SMatthew G. Knepley   }
158df032b82SMatthew G. Knepley   *gmsh_elems = elements;
159df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
160df032b82SMatthew G. Knepley }
1617d282ae0SMichael Lange 
162331830f3SMatthew G. Knepley /*@
1637d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
164331830f3SMatthew G. Knepley 
165331830f3SMatthew G. Knepley   Collective on comm
166331830f3SMatthew G. Knepley 
167331830f3SMatthew G. Knepley   Input Parameters:
168331830f3SMatthew G. Knepley + comm  - The MPI communicator
169331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
170331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
171331830f3SMatthew G. Knepley 
172331830f3SMatthew G. Knepley   Output Parameter:
173331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
174331830f3SMatthew G. Knepley 
175331830f3SMatthew G. Knepley   Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format
1763d51f2daSMichael Lange   and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format
177331830f3SMatthew G. Knepley 
178331830f3SMatthew G. Knepley   Level: beginner
179331830f3SMatthew G. Knepley 
180331830f3SMatthew G. Knepley .keywords: mesh,Gmsh
181331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
182331830f3SMatthew G. Knepley @*/
183331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
184331830f3SMatthew G. Knepley {
18504d1ad83SMichael Lange   PetscViewerType vtype;
1863c67d5adSSatish Balay   GmshElement   *gmsh_elem = NULL;
187331830f3SMatthew G. Knepley   PetscSection   coordSection;
188331830f3SMatthew G. Knepley   Vec            coordinates;
1896fbe17bfSStefano Zampini   PetscBT        periodicV = NULL, periodicC = NULL;
190331830f3SMatthew G. Knepley   PetscScalar   *coords, *coordsIn = NULL;
191f45c9589SStefano Zampini   PetscInt       dim = 0, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL;
192dc0ede3bSMatthew G. Knepley   int            i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum;
193331830f3SMatthew G. Knepley   PetscMPIInt    num_proc, rank;
194331830f3SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN];
195d3f73514SStefano Zampini   PetscBool      match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE;
196331830f3SMatthew G. Knepley   PetscErrorCode ierr;
197331830f3SMatthew G. Knepley 
198331830f3SMatthew G. Knepley   PetscFunctionBegin;
199331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
200331830f3SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr);
201331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
202331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
2033b3bc66dSMichael Lange   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
20404d1ad83SMichael Lange   ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr);
20504d1ad83SMichael Lange   ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
206d3f73514SStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_periodic",&periodic,NULL);CHKERRQ(ierr);
207d3f73514SStefano Zampini   ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_usemarker",&usemarker,NULL);CHKERRQ(ierr);
208abc86ac4SMichael Lange   if (!rank || binary) {
209331830f3SMatthew G. Knepley     PetscBool match;
210331830f3SMatthew G. Knepley     int       fileType, dataSize;
211f6ac7a6aSMichael Lange     float     version;
212331830f3SMatthew G. Knepley 
213331830f3SMatthew G. Knepley     /* Read header */
214060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
21504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
216331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
217060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
218f6ac7a6aSMichael Lange     snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
219f6ac7a6aSMichael Lange     if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line);
220f6ac7a6aSMichael Lange     if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0");
221331830f3SMatthew 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);
22204d1ad83SMichael Lange     if (binary) {
223b9eae255SMichael Lange       int checkInt;
224060da220SMatthew G. Knepley       ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
22504d1ad83SMichael Lange       if (checkInt != 1) {
226b9eae255SMichael Lange         ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr);
22704d1ad83SMichael Lange         if (checkInt == 1) bswap = PETSC_TRUE;
22804d1ad83SMichael Lange         else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType);
22904d1ad83SMichael Lange       }
2300877b519SMichael Lange     } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType);
231060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
23204d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
233331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
234dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
235060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
236dc0ede3bSMatthew G. Knepley     ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
237dc0ede3bSMatthew G. Knepley     if (match) {
238dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
239dc0ede3bSMatthew G. Knepley       snum = sscanf(line, "%d", &numRegions);
240dc0ede3bSMatthew G. Knepley       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
241a8667449SSander Arens       for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);}
242dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
243dc0ede3bSMatthew G. Knepley       ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
244dc0ede3bSMatthew G. Knepley       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
245dc0ede3bSMatthew G. Knepley       /* Initial read for vertex section */
246dc0ede3bSMatthew G. Knepley       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
247dc0ede3bSMatthew G. Knepley     }
248dc0ede3bSMatthew G. Knepley     /* Read vertices */
24904d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
250331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
251060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
2520877b519SMichael Lange     snum = sscanf(line, "%d", &numVertices);
2530877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
254854ce69bSBarry Smith     ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr);
25513aecfe5SMichael Lange     if (binary) {
25613aecfe5SMichael Lange       size_t doubleSize, intSize;
25713aecfe5SMichael Lange       PetscInt elementSize;
25813aecfe5SMichael Lange       char *buffer;
25913aecfe5SMichael Lange       PetscScalar *baseptr;
26013aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr);
26113aecfe5SMichael Lange       ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr);
26213aecfe5SMichael Lange       elementSize = (intSize + 3*doubleSize);
26313aecfe5SMichael Lange       ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr);
26413aecfe5SMichael Lange       ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr);
26513aecfe5SMichael Lange       if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr);
266331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
26713aecfe5SMichael Lange         baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize));
26813aecfe5SMichael Lange         coordsIn[v*3+0] = baseptr[0];
26913aecfe5SMichael Lange         coordsIn[v*3+1] = baseptr[1];
27013aecfe5SMichael Lange         coordsIn[v*3+2] = baseptr[2];
27113aecfe5SMichael Lange       }
27213aecfe5SMichael Lange       ierr = PetscFree(buffer);CHKERRQ(ierr);
27313aecfe5SMichael Lange     } else {
27413aecfe5SMichael Lange       for (v = 0; v < numVertices; ++v) {
275060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
276060da220SMatthew G. Knepley         ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
277b9eae255SMichael Lange         if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1);
278331830f3SMatthew G. Knepley       }
27913aecfe5SMichael Lange     }
280060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28104d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
282331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
283331830f3SMatthew G. Knepley     /* Read cells */
284060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
28504d1ad83SMichael Lange     ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr);
286331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
287060da220SMatthew G. Knepley     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);;
2880877b519SMichael Lange     snum = sscanf(line, "%d", &numCells);
2890877b519SMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
290331830f3SMatthew G. Knepley   }
291db397164SMichael Lange 
292abc86ac4SMichael Lange   if (!rank || binary) {
293a4bb7517SMichael Lange     /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
294a4bb7517SMichael Lange        file contents multiple times to figure out the true number of cells and facets
295a4bb7517SMichael Lange        in the given mesh. To make this more efficient we read the file contents only
296a4bb7517SMichael Lange        once and store them in memory, while determining the true number of cells. */
2973d51f2daSMichael Lange     ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr);
298a4bb7517SMichael Lange     for (trueNumCells=0, c = 0; c < numCells; ++c) {
299a4bb7517SMichael Lange       if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;}
300a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) trueNumCells++;
301db397164SMichael Lange     }
302d08df55aSStefano Zampini     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
3031b82c3ebSMatthew G. Knepley     ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr);
304331830f3SMatthew G. Knepley     if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
305d08df55aSStefano Zampini     if (periodic) {
306f45c9589SStefano Zampini       PetscInt pVert, *periodicMapT, *aux;
307d08df55aSStefano Zampini       int      numPeriodic;
308d08df55aSStefano Zampini 
309d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
310d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr);
311d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
312f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr);
3136fbe17bfSStefano Zampini       ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr);
314f45c9589SStefano Zampini       for (i = 0; i < numVertices; i++) periodicMapT[i] = i;
315d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
316d08df55aSStefano Zampini       snum = sscanf(line, "%d", &numPeriodic);
317d08df55aSStefano Zampini       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
318d08df55aSStefano Zampini       for (i = 0; i < numPeriodic; i++) {
319d08df55aSStefano Zampini         int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes;
320d08df55aSStefano Zampini 
321d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
322d08df55aSStefano Zampini         snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */
323d08df55aSStefano Zampini         if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
32472ffbcc9SStefano Zampini         /* Type of tranformation, discarded */
32572ffbcc9SStefano Zampini         ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
326d08df55aSStefano Zampini         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
327d08df55aSStefano Zampini         snum = sscanf(line, "%d", &nNodes);
328d08df55aSStefano Zampini         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
329d08df55aSStefano Zampini         for (j = 0; j < nNodes; j++) {
330d08df55aSStefano Zampini           PetscInt slaveNode, masterNode;
331d08df55aSStefano Zampini 
332d08df55aSStefano Zampini           ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
333d08df55aSStefano Zampini           snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
334d08df55aSStefano Zampini           if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
335f45c9589SStefano Zampini           periodicMapT[slaveNode - 1] = masterNode - 1;
336437bdd5bSStefano Zampini           ierr = PetscBTSet(periodicV, slaveNode - 1);CHKERRQ(ierr);
337437bdd5bSStefano Zampini           ierr = PetscBTSet(periodicV, masterNode - 1);CHKERRQ(ierr);
338d08df55aSStefano Zampini         }
339d08df55aSStefano Zampini       }
340d08df55aSStefano Zampini       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
341d08df55aSStefano Zampini       ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr);
342d08df55aSStefano Zampini       if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file");
343d08df55aSStefano Zampini       /* we may have slaves of slaves */
344d08df55aSStefano Zampini       for (i = 0; i < numVertices; i++) {
345f45c9589SStefano Zampini         while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) {
346f45c9589SStefano Zampini           periodicMapT[i] = periodicMapT[periodicMapT[i]];
347d08df55aSStefano Zampini         }
348d08df55aSStefano Zampini       }
349f45c9589SStefano Zampini       /* periodicMap : from old to new numbering (periodic vertices excluded)
350f45c9589SStefano Zampini          periodicMapI: from new to old numbering */
351f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr);
352f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr);
353f45c9589SStefano Zampini       ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr);
354f45c9589SStefano Zampini       for (i = 0, pVert = 0; i < numVertices; i++) {
355f45c9589SStefano Zampini         if (periodicMapT[i] != i) {
356f45c9589SStefano Zampini           pVert++;
357f45c9589SStefano Zampini         } else {
358f45c9589SStefano Zampini           aux[i] = i - pVert;
359f45c9589SStefano Zampini           periodicMapI[i - pVert] = i;
360f45c9589SStefano Zampini         }
361f45c9589SStefano Zampini       }
362f45c9589SStefano Zampini       for (i = 0 ; i < numVertices; i++) {
363f45c9589SStefano Zampini         periodicMap[i] = aux[periodicMapT[i]];
364f45c9589SStefano Zampini       }
365f45c9589SStefano Zampini       ierr = PetscFree(periodicMapT);CHKERRQ(ierr);
366f45c9589SStefano Zampini       ierr = PetscFree(aux);CHKERRQ(ierr);
367f45c9589SStefano Zampini       /* remove periodic vertices */
368f45c9589SStefano Zampini       numVertices = numVertices - pVert;
369d08df55aSStefano Zampini     }
370331830f3SMatthew G. Knepley   }
371abc86ac4SMichael Lange   /* For binary we read on all ranks, but only build the plex on rank 0 */
372abc86ac4SMichael Lange   if (binary && rank) {trueNumCells = 0; numVertices = 0;};
373a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
374331830f3SMatthew G. Knepley   ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr);
375331830f3SMatthew G. Knepley   if (!rank) {
376a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
377a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
378a4bb7517SMichael Lange         ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr);
379a4bb7517SMichael Lange         cell++;
380331830f3SMatthew G. Knepley       }
381331830f3SMatthew G. Knepley     }
382331830f3SMatthew G. Knepley   }
383331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
384a4bb7517SMichael Lange   /* Add cell-vertex connections */
385331830f3SMatthew G. Knepley   if (!rank) {
386331830f3SMatthew G. Knepley     PetscInt pcone[8], corner;
387437bdd5bSStefano Zampini 
388a4bb7517SMichael Lange     for (cell = 0, c = 0; c < numCells; ++c) {
389a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim) {
390a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
3916fbe17bfSStefano Zampini           const PetscInt cc = gmsh_elem[c].nodes[corner] - 1;
3926fbe17bfSStefano Zampini           pcone[corner] = periodicMap ? periodicMap[cc] + trueNumCells
3936fbe17bfSStefano Zampini                                       : cc + trueNumCells;
394331830f3SMatthew G. Knepley         }
39597ed6be6Ssarens         if (dim == 3) {
39697ed6be6Ssarens           /* Tetrahedra are inverted */
39797ed6be6Ssarens           if (gmsh_elem[c].numNodes == 4) {
39897ed6be6Ssarens             PetscInt tmp = pcone[0];
39997ed6be6Ssarens             pcone[0] = pcone[1];
40097ed6be6Ssarens             pcone[1] = tmp;
40197ed6be6Ssarens           }
40297ed6be6Ssarens           /* Hexahedra are inverted */
40397ed6be6Ssarens           if (gmsh_elem[c].numNodes == 8) {
40497ed6be6Ssarens             PetscInt tmp = pcone[1];
40597ed6be6Ssarens             pcone[1] = pcone[3];
40697ed6be6Ssarens             pcone[3] = tmp;
40797ed6be6Ssarens           }
40897ed6be6Ssarens         }
409a4bb7517SMichael Lange         ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr);
410a4bb7517SMichael Lange         cell++;
411331830f3SMatthew G. Knepley       }
412a4bb7517SMichael Lange     }
413331830f3SMatthew G. Knepley   }
4146228fc9fSJed Brown   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
415c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
416331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
417331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
418331830f3SMatthew G. Knepley   if (interpolate) {
419e56d480eSMatthew G. Knepley     DM idm = NULL;
420331830f3SMatthew G. Knepley 
421331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
422331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
423331830f3SMatthew G. Knepley     *dm  = idm;
424331830f3SMatthew G. Knepley   }
425d3f73514SStefano Zampini   if (usemarker) {
426d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
427d3f73514SStefano Zampini 
428d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm,"marker");CHKERRQ(ierr);
429d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
430d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
431d3f73514SStefano Zampini       PetscInt suppSize;
432d3f73514SStefano Zampini 
433d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
434d3f73514SStefano Zampini       if (suppSize == 1) {
435d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
436d3f73514SStefano Zampini 
437d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
438d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
439d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
440d3f73514SStefano Zampini         }
441d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
442d3f73514SStefano Zampini       }
443d3f73514SStefano Zampini     }
444d3f73514SStefano Zampini   }
44516ad7e67SMichael Lange 
44616ad7e67SMichael Lange   if (!rank) {
44716ad7e67SMichael Lange     /* Apply boundary IDs by finding the relevant facets with vertex joins */
44816ad7e67SMichael Lange     PetscInt pcone[8], corner, vStart, vEnd;
449d1a54cefSMatthew G. Knepley 
45016ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
45116ad7e67SMichael Lange     for (c = 0; c < numCells; ++c) {
452a4bb7517SMichael Lange       if (gmsh_elem[c].dim == dim-1) {
45316ad7e67SMichael Lange         const PetscInt *join;
454437bdd5bSStefano Zampini         PetscInt        joinSize;
455437bdd5bSStefano Zampini 
456a4bb7517SMichael Lange         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
457437bdd5bSStefano Zampini           const PetscInt cc = gmsh_elem[c].nodes[corner] - 1;
458437bdd5bSStefano Zampini           pcone[corner] = periodicMap ? periodicMap[cc] + vStart
459437bdd5bSStefano Zampini                                       : cc + vStart;
46016ad7e67SMichael Lange         }
461a4bb7517SMichael Lange         ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
462a4bb7517SMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id);
463c58f1c22SToby Isaac         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr);
464a4bb7517SMichael Lange         ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr);
46516ad7e67SMichael Lange       }
466a4bb7517SMichael Lange     }
4676e1dd89aSLawrence Mitchell 
4686e1dd89aSLawrence Mitchell     /* Create cell sets */
4696e1dd89aSLawrence Mitchell     for (cell = 0, c = 0; c < numCells; ++c) {
4706e1dd89aSLawrence Mitchell       if (gmsh_elem[c].dim == dim) {
4716e1dd89aSLawrence Mitchell         if (gmsh_elem[c].numTags > 0) {
4726e1dd89aSLawrence Mitchell           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4736e1dd89aSLawrence Mitchell           cell++;
4746e1dd89aSLawrence Mitchell         }
4756e1dd89aSLawrence Mitchell       }
4766e1dd89aSLawrence Mitchell     }
4770c070f12SSander Arens 
4780c070f12SSander Arens     /* Create vertex sets */
4790c070f12SSander Arens     for (c = 0; c < numCells; ++c) {
4800c070f12SSander Arens       if (gmsh_elem[c].dim == 0) {
4810c070f12SSander Arens         if (gmsh_elem[c].numTags > 0) {
482d08df55aSStefano Zampini           PetscInt vid =  periodicMap ? periodicMap[gmsh_elem[c].nodes[0] - 1] + vStart
483d08df55aSStefano Zampini                                       : gmsh_elem[c].nodes[0] - 1 + vStart;
484d08df55aSStefano Zampini           ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr);
4850c070f12SSander Arens         }
4860c070f12SSander Arens       }
4870c070f12SSander Arens     }
48816ad7e67SMichael Lange   }
48916ad7e67SMichael Lange 
490331830f3SMatthew G. Knepley   /* Read coordinates */
491979c4b60SMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
492331830f3SMatthew G. Knepley   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
493331830f3SMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
494f45c9589SStefano Zampini   if (periodicMap) { /* we need to localize coordinates on cells */
495f45c9589SStefano Zampini     ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr);
496f45c9589SStefano Zampini   } else {
49775b5763bSMichael Lange     ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr);
498f45c9589SStefano Zampini   }
49975b5763bSMichael Lange   for (v = trueNumCells; v < trueNumCells+numVertices; ++v) {
500331830f3SMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
501331830f3SMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
502331830f3SMatthew G. Knepley   }
503f45c9589SStefano Zampini   if (periodicMap) {
504437bdd5bSStefano Zampini     ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr);
505f45c9589SStefano Zampini     for (cell = 0, c = 0; c < numCells; ++c) {
506f45c9589SStefano Zampini       if (gmsh_elem[c].dim == dim) {
507437bdd5bSStefano Zampini         PetscBool pc = PETSC_FALSE;
508437bdd5bSStefano Zampini         PetscInt  corner;
509437bdd5bSStefano Zampini         for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
510*3548e9b6SStefano Zampini           pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - 1));
511437bdd5bSStefano Zampini         }
512437bdd5bSStefano Zampini         if (pc) {
513437bdd5bSStefano Zampini           ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr);
514f45c9589SStefano Zampini           ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
515f45c9589SStefano Zampini           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr);
5166fbe17bfSStefano Zampini         }
517f45c9589SStefano Zampini         cell++;
518f45c9589SStefano Zampini       }
519f45c9589SStefano Zampini     }
520f45c9589SStefano Zampini   }
521331830f3SMatthew G. Knepley   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
522331830f3SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
5238b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
524331830f3SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
525331830f3SMatthew G. Knepley   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
5267635fff5Ssarens   ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr);
527331830f3SMatthew G. Knepley   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
528331830f3SMatthew G. Knepley   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
529331830f3SMatthew G. Knepley   if (!rank) {
530f45c9589SStefano Zampini 
531f45c9589SStefano Zampini     if (periodicMap) {
532f45c9589SStefano Zampini       PetscInt off;
533f45c9589SStefano Zampini 
534f45c9589SStefano Zampini       for (cell = 0, c = 0; c < numCells; ++c) {
535f45c9589SStefano Zampini         PetscInt pcone[8], corner;
536f45c9589SStefano Zampini         if (gmsh_elem[c].dim == dim) {
5376fbe17bfSStefano Zampini           if (PetscUnlikely(PetscBTLookup(periodicC,cell))) {
538f45c9589SStefano Zampini             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
539f45c9589SStefano Zampini               pcone[corner] = gmsh_elem[c].nodes[corner] - 1;
540f45c9589SStefano Zampini             }
541f45c9589SStefano Zampini             if (dim == 3) {
542f45c9589SStefano Zampini               /* Tetrahedra are inverted */
543f45c9589SStefano Zampini               if (gmsh_elem[c].numNodes == 4) {
544f45c9589SStefano Zampini                 PetscInt tmp = pcone[0];
545f45c9589SStefano Zampini                 pcone[0] = pcone[1];
546f45c9589SStefano Zampini                 pcone[1] = tmp;
547f45c9589SStefano Zampini               }
548f45c9589SStefano Zampini               /* Hexahedra are inverted */
549f45c9589SStefano Zampini               if (gmsh_elem[c].numNodes == 8) {
550f45c9589SStefano Zampini                 PetscInt tmp = pcone[1];
551f45c9589SStefano Zampini                 pcone[1] = pcone[3];
552f45c9589SStefano Zampini                 pcone[3] = tmp;
553f45c9589SStefano Zampini               }
554f45c9589SStefano Zampini             }
555f45c9589SStefano Zampini             ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
556f45c9589SStefano Zampini             for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) {
557f45c9589SStefano Zampini               PetscInt v = pcone[corner];
558f45c9589SStefano Zampini               for (d = 0; d < dim; ++d) {
559f45c9589SStefano Zampini                 coords[off++] = coordsIn[v*3+d];
560f45c9589SStefano Zampini               }
561f45c9589SStefano Zampini             }
5626fbe17bfSStefano Zampini           }
563f45c9589SStefano Zampini           cell++;
564f45c9589SStefano Zampini         }
565f45c9589SStefano Zampini       }
566f45c9589SStefano Zampini       for (v = 0; v < numVertices; ++v) {
567f45c9589SStefano Zampini         ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr);
568f45c9589SStefano Zampini         for (d = 0; d < dim; ++d) {
569f45c9589SStefano Zampini           coords[off+d] = coordsIn[periodicMapI[v]*3+d];
570f45c9589SStefano Zampini         }
571f45c9589SStefano Zampini       }
572f45c9589SStefano Zampini     } else {
573331830f3SMatthew G. Knepley       for (v = 0; v < numVertices; ++v) {
574331830f3SMatthew G. Knepley         for (d = 0; d < dim; ++d) {
575331830f3SMatthew G. Knepley           coords[v*dim+d] = coordsIn[v*3+d];
576331830f3SMatthew G. Knepley         }
577331830f3SMatthew G. Knepley       }
578331830f3SMatthew G. Knepley     }
579f45c9589SStefano Zampini   }
58090b157c4SStefano Zampini   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
581331830f3SMatthew G. Knepley   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
582331830f3SMatthew G. Knepley   ierr = PetscFree(coordsIn);CHKERRQ(ierr);
583331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
584331830f3SMatthew G. Knepley   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
585d08df55aSStefano Zampini   ierr = PetscFree(periodicMap);CHKERRQ(ierr);
586fcd9ca0aSStefano Zampini   ierr = PetscFree(periodicMapI);CHKERRQ(ierr);
5876fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr);
5886fbe17bfSStefano Zampini   ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr);
589a4bb7517SMichael Lange   /* Clean up intermediate storage */
59029403c87SMichael Lange   if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr);
5913b3bc66dSMichael Lange   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr);
592331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
593331830f3SMatthew G. Knepley }
594