xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 066ea43f7f75752f012be6cd06b6107ebe84cc6d)
1331830f3SMatthew G. Knepley #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3730356f6SLisandro Dalcin #include <petsc/private/hashmapi.h>
4331830f3SMatthew G. Knepley 
5*066ea43fSLisandro Dalcin #include <../src/dm/impls/plex/gmshlex.h>
6*066ea43fSLisandro Dalcin 
7*066ea43fSLisandro Dalcin #define GMSH_LEXORDER_ITEM(T, p)                                   \
8*066ea43fSLisandro Dalcin static int *GmshLexOrder_##T##_##p(void)                           \
9*066ea43fSLisandro Dalcin {                                                                  \
10*066ea43fSLisandro Dalcin   static int Gmsh_LexOrder_##T##_##p[GmshNumNodes_##T(p)] = {-1};  \
11*066ea43fSLisandro Dalcin   int *lex = Gmsh_LexOrder_##T##_##p;                              \
12*066ea43fSLisandro Dalcin   if (lex[0] == -1) (void)GmshLexOrder_##T(p, lex, 0);             \
13*066ea43fSLisandro Dalcin   return lex;                                                      \
14*066ea43fSLisandro Dalcin }
15*066ea43fSLisandro Dalcin 
16*066ea43fSLisandro Dalcin #define GMSH_LEXORDER_LIST(T) \
17*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  1)     \
18*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  2)     \
19*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  3)     \
20*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  4)     \
21*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  5)     \
22*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  6)     \
23*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  7)     \
24*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  8)     \
25*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  9)     \
26*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 10)
27*066ea43fSLisandro Dalcin 
28*066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(VTX, 0)
29*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(SEG)
30*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TRI)
31*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(QUA)
32*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TET)
33*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(HEX)
34*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PRI)
35*066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PYR)
36*066ea43fSLisandro Dalcin 
37*066ea43fSLisandro Dalcin typedef enum {
38*066ea43fSLisandro Dalcin   GMSH_VTX = 0,
39*066ea43fSLisandro Dalcin   GMSH_SEG = 1,
40*066ea43fSLisandro Dalcin   GMSH_TRI = 2,
41*066ea43fSLisandro Dalcin   GMSH_QUA = 3,
42*066ea43fSLisandro Dalcin   GMSH_TET = 4,
43*066ea43fSLisandro Dalcin   GMSH_HEX = 5,
44*066ea43fSLisandro Dalcin   GMSH_PRI = 6,
45*066ea43fSLisandro Dalcin   GMSH_PYR = 7,
46*066ea43fSLisandro Dalcin   GMSH_NUM_POLYTOPES = 8
47*066ea43fSLisandro Dalcin } GmshPolytopeType;
48*066ea43fSLisandro Dalcin 
49de78e4feSLisandro Dalcin typedef struct {
500598e1a2SLisandro Dalcin   int   cellType;
51*066ea43fSLisandro Dalcin   int   polytope;
520598e1a2SLisandro Dalcin   int   dim;
530598e1a2SLisandro Dalcin   int   order;
54*066ea43fSLisandro Dalcin   int   numVerts;
550598e1a2SLisandro Dalcin   int   numNodes;
56*066ea43fSLisandro Dalcin   int* (*lexorder)(void);
570598e1a2SLisandro Dalcin } GmshCellInfo;
580598e1a2SLisandro Dalcin 
59*066ea43fSLisandro Dalcin #define GmshCellEntry(cellType, polytope, dim, order) \
60*066ea43fSLisandro Dalcin   {cellType, GMSH_##polytope, dim, order, \
61*066ea43fSLisandro Dalcin    GmshNumNodes_##polytope(1), \
62*066ea43fSLisandro Dalcin    GmshNumNodes_##polytope(order), \
63*066ea43fSLisandro Dalcin    GmshLexOrder_##polytope##_##order}
640598e1a2SLisandro Dalcin 
650598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = {
66*066ea43fSLisandro Dalcin   GmshCellEntry(  15, VTX, 0,  0),
670598e1a2SLisandro Dalcin 
68*066ea43fSLisandro Dalcin   GmshCellEntry(   1, SEG, 1,  1),
69*066ea43fSLisandro Dalcin   GmshCellEntry(   8, SEG, 1,  2),
70*066ea43fSLisandro Dalcin   GmshCellEntry(  26, SEG, 1,  3),
71*066ea43fSLisandro Dalcin   GmshCellEntry(  27, SEG, 1,  4),
72*066ea43fSLisandro Dalcin   GmshCellEntry(  28, SEG, 1,  5),
73*066ea43fSLisandro Dalcin   GmshCellEntry(  62, SEG, 1,  6),
74*066ea43fSLisandro Dalcin   GmshCellEntry(  63, SEG, 1,  7),
75*066ea43fSLisandro Dalcin   GmshCellEntry(  64, SEG, 1,  8),
76*066ea43fSLisandro Dalcin   GmshCellEntry(  65, SEG, 1,  9),
77*066ea43fSLisandro Dalcin   GmshCellEntry(  66, SEG, 1, 10),
780598e1a2SLisandro Dalcin 
79*066ea43fSLisandro Dalcin   GmshCellEntry(   2, TRI, 2,  1),
80*066ea43fSLisandro Dalcin   GmshCellEntry(   9, TRI, 2,  2),
81*066ea43fSLisandro Dalcin   GmshCellEntry(  21, TRI, 2,  3),
82*066ea43fSLisandro Dalcin   GmshCellEntry(  23, TRI, 2,  4),
83*066ea43fSLisandro Dalcin   GmshCellEntry(  25, TRI, 2,  5),
84*066ea43fSLisandro Dalcin   GmshCellEntry(  42, TRI, 2,  6),
85*066ea43fSLisandro Dalcin   GmshCellEntry(  43, TRI, 2,  7),
86*066ea43fSLisandro Dalcin   GmshCellEntry(  44, TRI, 2,  8),
87*066ea43fSLisandro Dalcin   GmshCellEntry(  45, TRI, 2,  9),
88*066ea43fSLisandro Dalcin   GmshCellEntry(  46, TRI, 2, 10),
890598e1a2SLisandro Dalcin 
90*066ea43fSLisandro Dalcin   GmshCellEntry(   3, QUA, 2,  1),
91*066ea43fSLisandro Dalcin   GmshCellEntry(  10, QUA, 2,  2),
92*066ea43fSLisandro Dalcin   GmshCellEntry(  36, QUA, 2,  3),
93*066ea43fSLisandro Dalcin   GmshCellEntry(  37, QUA, 2,  4),
94*066ea43fSLisandro Dalcin   GmshCellEntry(  38, QUA, 2,  5),
95*066ea43fSLisandro Dalcin   GmshCellEntry(  47, QUA, 2,  6),
96*066ea43fSLisandro Dalcin   GmshCellEntry(  48, QUA, 2,  7),
97*066ea43fSLisandro Dalcin   GmshCellEntry(  49, QUA, 2,  8),
98*066ea43fSLisandro Dalcin   GmshCellEntry(  50, QUA, 2,  9),
99*066ea43fSLisandro Dalcin   GmshCellEntry(  51, QUA, 2, 10),
1000598e1a2SLisandro Dalcin 
101*066ea43fSLisandro Dalcin   GmshCellEntry(   4, TET, 3,  1),
102*066ea43fSLisandro Dalcin   GmshCellEntry(  11, TET, 3,  2),
103*066ea43fSLisandro Dalcin   GmshCellEntry(  29, TET, 3,  3),
104*066ea43fSLisandro Dalcin   GmshCellEntry(  30, TET, 3,  4),
105*066ea43fSLisandro Dalcin   GmshCellEntry(  31, TET, 3,  5),
106*066ea43fSLisandro Dalcin   GmshCellEntry(  71, TET, 3,  6),
107*066ea43fSLisandro Dalcin   GmshCellEntry(  72, TET, 3,  7),
108*066ea43fSLisandro Dalcin   GmshCellEntry(  73, TET, 3,  8),
109*066ea43fSLisandro Dalcin   GmshCellEntry(  74, TET, 3,  9),
110*066ea43fSLisandro Dalcin   GmshCellEntry(  75, TET, 3, 10),
1110598e1a2SLisandro Dalcin 
112*066ea43fSLisandro Dalcin   GmshCellEntry(   5, HEX, 3,  1),
113*066ea43fSLisandro Dalcin   GmshCellEntry(  12, HEX, 3,  2),
114*066ea43fSLisandro Dalcin   GmshCellEntry(  92, HEX, 3,  3),
115*066ea43fSLisandro Dalcin   GmshCellEntry(  93, HEX, 3,  4),
116*066ea43fSLisandro Dalcin   GmshCellEntry(  94, HEX, 3,  5),
117*066ea43fSLisandro Dalcin   GmshCellEntry(  95, HEX, 3,  6),
118*066ea43fSLisandro Dalcin   GmshCellEntry(  96, HEX, 3,  7),
119*066ea43fSLisandro Dalcin   GmshCellEntry(  97, HEX, 3,  8),
120*066ea43fSLisandro Dalcin   GmshCellEntry(  98, HEX, 3,  9),
121*066ea43fSLisandro Dalcin   GmshCellEntry(  -1, HEX, 3, 10),
1220598e1a2SLisandro Dalcin 
123*066ea43fSLisandro Dalcin   GmshCellEntry(   6, PRI, 3,  1),
124*066ea43fSLisandro Dalcin   GmshCellEntry(  13, PRI, 3,  2),
125*066ea43fSLisandro Dalcin   GmshCellEntry(  90, PRI, 3,  3),
126*066ea43fSLisandro Dalcin   GmshCellEntry(  91, PRI, 3,  4),
127*066ea43fSLisandro Dalcin   GmshCellEntry( 106, PRI, 3,  5),
128*066ea43fSLisandro Dalcin   GmshCellEntry( 107, PRI, 3,  6),
129*066ea43fSLisandro Dalcin   GmshCellEntry( 108, PRI, 3,  7),
130*066ea43fSLisandro Dalcin   GmshCellEntry( 109, PRI, 3,  8),
131*066ea43fSLisandro Dalcin   GmshCellEntry( 110, PRI, 3,  9),
132*066ea43fSLisandro Dalcin   GmshCellEntry(  -1, PRI, 3, 10),
1330598e1a2SLisandro Dalcin 
134*066ea43fSLisandro Dalcin   GmshCellEntry(   7, PYR, 3,  1),
135*066ea43fSLisandro Dalcin   GmshCellEntry(  14, PYR, 3,  2),
136*066ea43fSLisandro Dalcin   GmshCellEntry( 118, PYR, 3,  3),
137*066ea43fSLisandro Dalcin   GmshCellEntry( 119, PYR, 3,  4),
138*066ea43fSLisandro Dalcin   GmshCellEntry( 120, PYR, 3,  5),
139*066ea43fSLisandro Dalcin   GmshCellEntry( 121, PYR, 3,  6),
140*066ea43fSLisandro Dalcin   GmshCellEntry( 122, PYR, 3,  7),
141*066ea43fSLisandro Dalcin   GmshCellEntry( 123, PYR, 3,  8),
142*066ea43fSLisandro Dalcin   GmshCellEntry( 124, PYR, 3,  9),
143*066ea43fSLisandro Dalcin   GmshCellEntry(  -1, PYR, 3, 10)
1440598e1a2SLisandro Dalcin 
1450598e1a2SLisandro Dalcin #if 0
146*066ea43fSLisandro Dalcin   { 20, GMSH_TRI, 2, 3, 3,  9, NULL},
147*066ea43fSLisandro Dalcin   { 16, GMSH_QUA, 2, 2, 4,  8, NULL},
148*066ea43fSLisandro Dalcin   { 17, GMSH_HEX, 3, 2, 8, 20, NULL},
149*066ea43fSLisandro Dalcin   { 18, GMSH_PRI, 3, 2, 6, 15, NULL},
150*066ea43fSLisandro Dalcin   { 19, GMSH_PYR, 3, 2, 5, 13, NULL},
1510598e1a2SLisandro Dalcin #endif
1520598e1a2SLisandro Dalcin };
1530598e1a2SLisandro Dalcin 
1540598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150];
1550598e1a2SLisandro Dalcin 
1560598e1a2SLisandro Dalcin static PetscErrorCode GmshCellInfoSetUp(void)
1570598e1a2SLisandro Dalcin {
1580598e1a2SLisandro Dalcin   size_t           i, n;
1590598e1a2SLisandro Dalcin   static PetscBool called = PETSC_FALSE;
1600598e1a2SLisandro Dalcin 
1610598e1a2SLisandro Dalcin   if (called) return 0;
1620598e1a2SLisandro Dalcin   PetscFunctionBegin;
1630598e1a2SLisandro Dalcin   called = PETSC_TRUE;
1640598e1a2SLisandro Dalcin   n = sizeof(GmshCellMap)/sizeof(GmshCellMap[0]);
1650598e1a2SLisandro Dalcin   for (i = 0; i < n; ++i) {
1660598e1a2SLisandro Dalcin     GmshCellMap[i].cellType = -1;
167*066ea43fSLisandro Dalcin     GmshCellMap[i].polytope = -1;
1680598e1a2SLisandro Dalcin   }
1690598e1a2SLisandro Dalcin   n = sizeof(GmshCellTable)/sizeof(GmshCellTable[0]);
170*066ea43fSLisandro Dalcin   for (i = 0; i < n; ++i) {
171*066ea43fSLisandro Dalcin     if (GmshCellTable[i].cellType <= 0) continue;
172*066ea43fSLisandro Dalcin     GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i];
173*066ea43fSLisandro Dalcin   }
1740598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
1750598e1a2SLisandro Dalcin }
1760598e1a2SLisandro Dalcin 
1770598e1a2SLisandro Dalcin #define GmshCellTypeCheck(ct) 0; do { \
1780598e1a2SLisandro Dalcin     const int _ct_ = (int)ct; \
1790598e1a2SLisandro Dalcin     if (_ct_ < 0 || _ct_ >= (int)(sizeof(GmshCellMap)/sizeof(GmshCellMap[0]))) \
1800598e1a2SLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); \
1810598e1a2SLisandro Dalcin     if (GmshCellMap[_ct_].cellType != _ct_) \
1820598e1a2SLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \
183*066ea43fSLisandro Dalcin     if (GmshCellMap[_ct_].polytope == -1) \
1840598e1a2SLisandro Dalcin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \
1850598e1a2SLisandro Dalcin   } while(0)
1860598e1a2SLisandro Dalcin 
1870598e1a2SLisandro Dalcin 
1880598e1a2SLisandro Dalcin typedef struct {
189698a79b9SLisandro Dalcin   PetscViewer  viewer;
190698a79b9SLisandro Dalcin   int          fileFormat;
191698a79b9SLisandro Dalcin   int          dataSize;
192698a79b9SLisandro Dalcin   PetscBool    binary;
193698a79b9SLisandro Dalcin   PetscBool    byteSwap;
194698a79b9SLisandro Dalcin   size_t       wlen;
195698a79b9SLisandro Dalcin   void        *wbuf;
196698a79b9SLisandro Dalcin   size_t       slen;
197698a79b9SLisandro Dalcin   void        *sbuf;
1980598e1a2SLisandro Dalcin   PetscInt    *nbuf;
1990598e1a2SLisandro Dalcin   PetscInt     nodeStart;
2000598e1a2SLisandro Dalcin   PetscInt     nodeEnd;
20133a088b6SMatthew G. Knepley   PetscInt    *nodeMap;
202698a79b9SLisandro Dalcin } GmshFile;
203de78e4feSLisandro Dalcin 
204698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf)
205de78e4feSLisandro Dalcin {
206de78e4feSLisandro Dalcin   size_t         size = count * eltsize;
20711c56cb0SLisandro Dalcin   PetscErrorCode ierr;
20811c56cb0SLisandro Dalcin 
20911c56cb0SLisandro Dalcin   PetscFunctionBegin;
210698a79b9SLisandro Dalcin   if (gmsh->wlen < size) {
211698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr);
212698a79b9SLisandro Dalcin     ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr);
213698a79b9SLisandro Dalcin     gmsh->wlen = size;
214de78e4feSLisandro Dalcin   }
215698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->wbuf : NULL;
216de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
217de78e4feSLisandro Dalcin }
218de78e4feSLisandro Dalcin 
219698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf)
220698a79b9SLisandro Dalcin {
221698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
222698a79b9SLisandro Dalcin   size_t         size = count * dataSize;
223698a79b9SLisandro Dalcin   PetscErrorCode ierr;
224698a79b9SLisandro Dalcin 
225698a79b9SLisandro Dalcin   PetscFunctionBegin;
226698a79b9SLisandro Dalcin   if (gmsh->slen < size) {
227698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr);
228698a79b9SLisandro Dalcin     ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr);
229698a79b9SLisandro Dalcin     gmsh->slen = size;
230698a79b9SLisandro Dalcin   }
231698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->sbuf : NULL;
232698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
233698a79b9SLisandro Dalcin }
234698a79b9SLisandro Dalcin 
235698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype)
236de78e4feSLisandro Dalcin {
237de78e4feSLisandro Dalcin   PetscErrorCode ierr;
238de78e4feSLisandro Dalcin   PetscFunctionBegin;
239698a79b9SLisandro Dalcin   ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr);
240698a79b9SLisandro Dalcin   if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);}
241698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
242698a79b9SLisandro Dalcin }
243698a79b9SLisandro Dalcin 
244698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count)
245698a79b9SLisandro Dalcin {
246698a79b9SLisandro Dalcin   PetscErrorCode ierr;
247698a79b9SLisandro Dalcin   PetscFunctionBegin;
248698a79b9SLisandro Dalcin   ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr);
249698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
250698a79b9SLisandro Dalcin }
251698a79b9SLisandro Dalcin 
252d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match)
253d6689ca9SLisandro Dalcin {
254d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
255d6689ca9SLisandro Dalcin   PetscFunctionBegin;
256d6689ca9SLisandro Dalcin   ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr);
257d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
258d6689ca9SLisandro Dalcin }
259d6689ca9SLisandro Dalcin 
260d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN])
261d6689ca9SLisandro Dalcin {
262d6689ca9SLisandro Dalcin   PetscBool      match;
263d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
264d6689ca9SLisandro Dalcin 
265d6689ca9SLisandro Dalcin   PetscFunctionBegin;
266d6689ca9SLisandro Dalcin   ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr);
2670598e1a2SLisandro Dalcin   if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s",Section);
268d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
269d6689ca9SLisandro Dalcin }
270d6689ca9SLisandro Dalcin 
271d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN])
272d6689ca9SLisandro Dalcin {
273d6689ca9SLisandro Dalcin   PetscBool      match;
274d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
275d6689ca9SLisandro Dalcin 
276d6689ca9SLisandro Dalcin   PetscFunctionBegin;
277d6689ca9SLisandro Dalcin   while (PETSC_TRUE) {
278d6689ca9SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
279d6689ca9SLisandro Dalcin     ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr);
280d6689ca9SLisandro Dalcin     if (!match) break;
281d6689ca9SLisandro Dalcin     while (PETSC_TRUE) {
282d6689ca9SLisandro Dalcin       ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
283d6689ca9SLisandro Dalcin       ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr);
284d6689ca9SLisandro Dalcin       if (match) break;
285d6689ca9SLisandro Dalcin     }
286d6689ca9SLisandro Dalcin   }
287d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
288d6689ca9SLisandro Dalcin }
289d6689ca9SLisandro Dalcin 
290d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN])
291d6689ca9SLisandro Dalcin {
292d6689ca9SLisandro Dalcin   PetscErrorCode ierr;
293d6689ca9SLisandro Dalcin   PetscFunctionBegin;
294d6689ca9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
295d6689ca9SLisandro Dalcin   ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr);
296d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
297d6689ca9SLisandro Dalcin }
298d6689ca9SLisandro Dalcin 
299698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count)
300698a79b9SLisandro Dalcin {
301698a79b9SLisandro Dalcin   PetscInt       i;
302698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
303698a79b9SLisandro Dalcin   PetscErrorCode ierr;
304698a79b9SLisandro Dalcin 
305698a79b9SLisandro Dalcin   PetscFunctionBegin;
306698a79b9SLisandro Dalcin   if (dataSize == sizeof(PetscInt)) {
307698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr);
308698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(int)) {
309698a79b9SLisandro Dalcin     int *ibuf = NULL;
31089d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
311698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr);
312698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
313698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(long)) {
314698a79b9SLisandro Dalcin     long *ibuf = NULL;
31589d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
316698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr);
317698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
318698a79b9SLisandro Dalcin   } else if (dataSize == sizeof(PetscInt64)) {
319698a79b9SLisandro Dalcin     PetscInt64 *ibuf = NULL;
32089d5b078SSatish Balay     ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr);
321698a79b9SLisandro Dalcin     ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr);
322698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
323698a79b9SLisandro Dalcin   }
324698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
325698a79b9SLisandro Dalcin }
326698a79b9SLisandro Dalcin 
327698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count)
328698a79b9SLisandro Dalcin {
329698a79b9SLisandro Dalcin   PetscErrorCode ierr;
330698a79b9SLisandro Dalcin   PetscFunctionBegin;
331698a79b9SLisandro Dalcin   ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr);
332698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
333698a79b9SLisandro Dalcin }
334698a79b9SLisandro Dalcin 
335698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count)
336698a79b9SLisandro Dalcin {
337698a79b9SLisandro Dalcin   PetscErrorCode ierr;
338698a79b9SLisandro Dalcin   PetscFunctionBegin;
339698a79b9SLisandro Dalcin   ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr);
340de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
341de78e4feSLisandro Dalcin }
342de78e4feSLisandro Dalcin 
343de78e4feSLisandro Dalcin typedef struct {
3440598e1a2SLisandro Dalcin   PetscInt id;       /* Entity ID */
3450598e1a2SLisandro Dalcin   PetscInt dim;      /* Dimension */
346de78e4feSLisandro Dalcin   double   bbox[6];  /* Bounding box */
347de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
348de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
349de78e4feSLisandro Dalcin } GmshEntity;
350de78e4feSLisandro Dalcin 
351de78e4feSLisandro Dalcin typedef struct {
352730356f6SLisandro Dalcin   GmshEntity *entity[4];
353730356f6SLisandro Dalcin   PetscHMapI  entityMap[4];
354730356f6SLisandro Dalcin } GmshEntities;
355730356f6SLisandro Dalcin 
356698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities)
357730356f6SLisandro Dalcin {
358730356f6SLisandro Dalcin   PetscInt       dim;
359730356f6SLisandro Dalcin   PetscErrorCode ierr;
360730356f6SLisandro Dalcin 
361730356f6SLisandro Dalcin   PetscFunctionBegin;
362730356f6SLisandro Dalcin   ierr = PetscNew(entities);CHKERRQ(ierr);
363730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
364730356f6SLisandro Dalcin     ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr);
365730356f6SLisandro Dalcin     ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
366730356f6SLisandro Dalcin   }
367730356f6SLisandro Dalcin   PetscFunctionReturn(0);
368730356f6SLisandro Dalcin }
369730356f6SLisandro Dalcin 
3700598e1a2SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities)
3710598e1a2SLisandro Dalcin {
3720598e1a2SLisandro Dalcin   PetscInt       dim;
3730598e1a2SLisandro Dalcin   PetscErrorCode ierr;
3740598e1a2SLisandro Dalcin 
3750598e1a2SLisandro Dalcin   PetscFunctionBegin;
3760598e1a2SLisandro Dalcin   if (!*entities) PetscFunctionReturn(0);
3770598e1a2SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
3780598e1a2SLisandro Dalcin     ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr);
3790598e1a2SLisandro Dalcin     ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr);
3800598e1a2SLisandro Dalcin   }
3810598e1a2SLisandro Dalcin   ierr = PetscFree((*entities));CHKERRQ(ierr);
3820598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
3830598e1a2SLisandro Dalcin }
3840598e1a2SLisandro Dalcin 
385730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity)
386730356f6SLisandro Dalcin {
387730356f6SLisandro Dalcin   PetscErrorCode ierr;
3880598e1a2SLisandro Dalcin 
389730356f6SLisandro Dalcin   PetscFunctionBegin;
390730356f6SLisandro Dalcin   ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr);
391730356f6SLisandro Dalcin   entities->entity[dim][index].dim = dim;
392730356f6SLisandro Dalcin   entities->entity[dim][index].id  = eid;
393730356f6SLisandro Dalcin   if (entity) *entity = &entities->entity[dim][index];
394730356f6SLisandro Dalcin   PetscFunctionReturn(0);
395730356f6SLisandro Dalcin }
396730356f6SLisandro Dalcin 
397730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity)
398730356f6SLisandro Dalcin {
399730356f6SLisandro Dalcin   PetscInt       index;
400730356f6SLisandro Dalcin   PetscErrorCode ierr;
401730356f6SLisandro Dalcin 
402730356f6SLisandro Dalcin   PetscFunctionBegin;
403730356f6SLisandro Dalcin   ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr);
404730356f6SLisandro Dalcin   *entity = &entities->entity[dim][index];
405730356f6SLisandro Dalcin   PetscFunctionReturn(0);
406730356f6SLisandro Dalcin }
407730356f6SLisandro Dalcin 
4080598e1a2SLisandro Dalcin typedef struct {
4090598e1a2SLisandro Dalcin   PetscInt *id;   /* Node IDs */
4100598e1a2SLisandro Dalcin   double   *xyz;  /* Coordinates */
4110598e1a2SLisandro Dalcin } GmshNodes;
4120598e1a2SLisandro Dalcin 
4130598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes)
414730356f6SLisandro Dalcin {
415730356f6SLisandro Dalcin   PetscErrorCode ierr;
416730356f6SLisandro Dalcin 
417730356f6SLisandro Dalcin   PetscFunctionBegin;
4180598e1a2SLisandro Dalcin   ierr = PetscNew(nodes);CHKERRQ(ierr);
4190598e1a2SLisandro Dalcin   ierr = PetscMalloc1(count*1, &(*nodes)->id);CHKERRQ(ierr);
4200598e1a2SLisandro Dalcin   ierr = PetscMalloc1(count*3, &(*nodes)->xyz);CHKERRQ(ierr);
4210598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
422730356f6SLisandro Dalcin }
4230598e1a2SLisandro Dalcin 
4240598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes)
4250598e1a2SLisandro Dalcin {
4260598e1a2SLisandro Dalcin   PetscErrorCode ierr;
4270598e1a2SLisandro Dalcin   PetscFunctionBegin;
4280598e1a2SLisandro Dalcin   if (!*nodes) PetscFunctionReturn(0);
4290598e1a2SLisandro Dalcin   ierr = PetscFree((*nodes)->id);CHKERRQ(ierr);
4300598e1a2SLisandro Dalcin   ierr = PetscFree((*nodes)->xyz);CHKERRQ(ierr);
4310598e1a2SLisandro Dalcin   ierr = PetscFree((*nodes));CHKERRQ(ierr);
432730356f6SLisandro Dalcin   PetscFunctionReturn(0);
433730356f6SLisandro Dalcin }
434730356f6SLisandro Dalcin 
435730356f6SLisandro Dalcin typedef struct {
4360598e1a2SLisandro Dalcin   PetscInt id;       /* Element ID */
4370598e1a2SLisandro Dalcin   PetscInt dim;      /* Dimension */
438de78e4feSLisandro Dalcin   PetscInt cellType; /* Cell type */
4390598e1a2SLisandro Dalcin   PetscInt numVerts; /* Size of vertex array */
440de78e4feSLisandro Dalcin   PetscInt numNodes; /* Size of node array */
4410598e1a2SLisandro Dalcin   PetscInt *nodes;   /* Vertex/Node array */
442de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
443de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
444de78e4feSLisandro Dalcin } GmshElement;
445de78e4feSLisandro Dalcin 
4460598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements)
4470598e1a2SLisandro Dalcin {
4480598e1a2SLisandro Dalcin   PetscErrorCode ierr;
4490598e1a2SLisandro Dalcin 
4500598e1a2SLisandro Dalcin   PetscFunctionBegin;
4510598e1a2SLisandro Dalcin   ierr = PetscCalloc1(count, elements);CHKERRQ(ierr);
4520598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4530598e1a2SLisandro Dalcin }
4540598e1a2SLisandro Dalcin 
4550598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsDestroy(GmshElement **elements)
4560598e1a2SLisandro Dalcin {
4570598e1a2SLisandro Dalcin   PetscErrorCode ierr;
4580598e1a2SLisandro Dalcin 
4590598e1a2SLisandro Dalcin   PetscFunctionBegin;
4600598e1a2SLisandro Dalcin   if (!*elements) PetscFunctionReturn(0);
4610598e1a2SLisandro Dalcin   ierr = PetscFree(*elements);CHKERRQ(ierr);
4620598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4630598e1a2SLisandro Dalcin }
4640598e1a2SLisandro Dalcin 
4650598e1a2SLisandro Dalcin typedef struct {
4660598e1a2SLisandro Dalcin   PetscInt       dim;
467*066ea43fSLisandro Dalcin   PetscInt       order;
4680598e1a2SLisandro Dalcin   GmshEntities  *entities;
4690598e1a2SLisandro Dalcin   PetscInt       numNodes;
4700598e1a2SLisandro Dalcin   GmshNodes     *nodelist;
4710598e1a2SLisandro Dalcin   PetscInt       numElems;
4720598e1a2SLisandro Dalcin   GmshElement   *elements;
4730598e1a2SLisandro Dalcin   PetscInt       numVerts;
4740598e1a2SLisandro Dalcin   PetscInt       numCells;
4750598e1a2SLisandro Dalcin   PetscInt      *periodMap;
4760598e1a2SLisandro Dalcin   PetscInt      *vertexMap;
4770598e1a2SLisandro Dalcin   PetscSegBuffer segbuf;
4780598e1a2SLisandro Dalcin } GmshMesh;
4790598e1a2SLisandro Dalcin 
4800598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshCreate(GmshMesh **mesh)
4810598e1a2SLisandro Dalcin {
4820598e1a2SLisandro Dalcin   PetscErrorCode ierr;
4830598e1a2SLisandro Dalcin 
4840598e1a2SLisandro Dalcin   PetscFunctionBegin;
4850598e1a2SLisandro Dalcin   ierr = PetscNew(mesh);CHKERRQ(ierr);
4860598e1a2SLisandro Dalcin   ierr = PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf);CHKERRQ(ierr);
4870598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4880598e1a2SLisandro Dalcin }
4890598e1a2SLisandro Dalcin 
4900598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh)
4910598e1a2SLisandro Dalcin {
4920598e1a2SLisandro Dalcin   PetscErrorCode ierr;
4930598e1a2SLisandro Dalcin 
4940598e1a2SLisandro Dalcin   PetscFunctionBegin;
4950598e1a2SLisandro Dalcin   if (!*mesh) PetscFunctionReturn(0);
4960598e1a2SLisandro Dalcin   ierr = GmshEntitiesDestroy(&(*mesh)->entities);CHKERRQ(ierr);
4970598e1a2SLisandro Dalcin   ierr = GmshNodesDestroy(&(*mesh)->nodelist);CHKERRQ(ierr);
4980598e1a2SLisandro Dalcin   ierr = GmshElementsDestroy(&(*mesh)->elements);CHKERRQ(ierr);
4990598e1a2SLisandro Dalcin   ierr = PetscFree((*mesh)->periodMap);CHKERRQ(ierr);
5000598e1a2SLisandro Dalcin   ierr = PetscFree((*mesh)->vertexMap);CHKERRQ(ierr);
5010598e1a2SLisandro Dalcin   ierr = PetscSegBufferDestroy(&(*mesh)->segbuf);CHKERRQ(ierr);
5020598e1a2SLisandro Dalcin   ierr = PetscFree((*mesh));CHKERRQ(ierr);
5030598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
5040598e1a2SLisandro Dalcin }
5050598e1a2SLisandro Dalcin 
5060598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh)
507de78e4feSLisandro Dalcin {
508698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
509698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
510de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
5110598e1a2SLisandro Dalcin   int            n, num, nid, snum;
5120598e1a2SLisandro Dalcin   GmshNodes      *nodes;
513de78e4feSLisandro Dalcin   PetscErrorCode ierr;
514de78e4feSLisandro Dalcin 
515de78e4feSLisandro Dalcin   PetscFunctionBegin;
516de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
517698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
5180598e1a2SLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
5190598e1a2SLisandro Dalcin   ierr = GmshNodesCreate(num, &nodes);CHKERRQ(ierr);
5200598e1a2SLisandro Dalcin   mesh->numNodes = num;
5210598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
5220598e1a2SLisandro Dalcin   for (n = 0; n < num; ++n) {
5230598e1a2SLisandro Dalcin     double *xyz = nodes->xyz + n*3;
52411c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
52511c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
5260598e1a2SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
52711c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
5280598e1a2SLisandro Dalcin     nodes->id[n] = nid;
52911c56cb0SLisandro Dalcin   }
53011c56cb0SLisandro Dalcin   PetscFunctionReturn(0);
53111c56cb0SLisandro Dalcin }
53211c56cb0SLisandro Dalcin 
533de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
534de78e4feSLisandro Dalcin    file contents multiple times to figure out the true number of cells and facets
535de78e4feSLisandro Dalcin    in the given mesh. To make this more efficient we read the file contents only
536de78e4feSLisandro Dalcin    once and store them in memory, while determining the true number of cells. */
5370598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v22(GmshFile* gmsh, GmshMesh *mesh)
53811c56cb0SLisandro Dalcin {
539698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
540698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
541698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
542de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
5430598e1a2SLisandro Dalcin   int            i, c, p, num, ibuf[1+4+1000], snum;
5440598e1a2SLisandro Dalcin   int            cellType, numElem, numVerts, numNodes, numTags;
545df032b82SMatthew G. Knepley   GmshElement   *elements;
5460598e1a2SLisandro Dalcin   PetscInt      *nodeMap = gmsh->nodeMap;
547df032b82SMatthew G. Knepley   PetscErrorCode ierr;
548df032b82SMatthew G. Knepley 
549df032b82SMatthew G. Knepley   PetscFunctionBegin;
550feb237baSPierre Jolivet   ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
551698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
5520598e1a2SLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
5530598e1a2SLisandro Dalcin   ierr = GmshElementsCreate(num, &elements);CHKERRQ(ierr);
5540598e1a2SLisandro Dalcin   mesh->numElems = num;
5550598e1a2SLisandro Dalcin   mesh->elements = elements;
556698a79b9SLisandro Dalcin   for (c = 0; c < num;) {
55711c56cb0SLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
55811c56cb0SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
5590598e1a2SLisandro Dalcin 
5600598e1a2SLisandro Dalcin     cellType = binary ? ibuf[0] : ibuf[1];
5610598e1a2SLisandro Dalcin     numElem  = binary ? ibuf[1] : 1;
562df032b82SMatthew G. Knepley     numTags  = ibuf[2];
5630598e1a2SLisandro Dalcin 
5640598e1a2SLisandro Dalcin     ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr);
5650598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
5660598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
5670598e1a2SLisandro Dalcin 
568df032b82SMatthew G. Knepley     for (i = 0; i < numElem; ++i, ++c) {
5690598e1a2SLisandro Dalcin       GmshElement *element = elements + c;
5700598e1a2SLisandro Dalcin       const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off;
5710598e1a2SLisandro Dalcin       const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1;
5720598e1a2SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf+off, nint, NULL, PETSC_ENUM);CHKERRQ(ierr);
5730598e1a2SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf+off, PETSC_ENUM, nint);CHKERRQ(ierr);}
5740598e1a2SLisandro Dalcin       element->id  = id[0];
5750598e1a2SLisandro Dalcin       element->dim = GmshCellMap[cellType].dim;
5760598e1a2SLisandro Dalcin       element->cellType = cellType;
5770598e1a2SLisandro Dalcin       element->numVerts = numVerts;
5780598e1a2SLisandro Dalcin       element->numNodes = numNodes;
5790598e1a2SLisandro Dalcin       element->numTags  = PetscMin(numTags, 4);
5800598e1a2SLisandro Dalcin       ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr);
5810598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
5820598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = tags[p];
583df032b82SMatthew G. Knepley     }
584df032b82SMatthew G. Knepley   }
585df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
586df032b82SMatthew G. Knepley }
5877d282ae0SMichael Lange 
588de78e4feSLisandro Dalcin /*
589de78e4feSLisandro Dalcin $Entities
590de78e4feSLisandro Dalcin   numPoints(unsigned long) numCurves(unsigned long)
591de78e4feSLisandro Dalcin     numSurfaces(unsigned long) numVolumes(unsigned long)
592de78e4feSLisandro Dalcin   // points
593de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
594de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
595de78e4feSLisandro Dalcin     numPhysicals(unsigned long) phyisicalTag[...](int)
596de78e4feSLisandro Dalcin   ...
597de78e4feSLisandro Dalcin   // curves
598de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
599de78e4feSLisandro Dalcin      boxMaxX(double) boxMaxY(double) boxMaxZ(double)
600de78e4feSLisandro Dalcin      numPhysicals(unsigned long) physicalTag[...](int)
601de78e4feSLisandro Dalcin      numBREPVert(unsigned long) tagBREPVert[...](int)
602de78e4feSLisandro Dalcin   ...
603de78e4feSLisandro Dalcin   // surfaces
604de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
605de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
606de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
607de78e4feSLisandro Dalcin     numBREPCurve(unsigned long) tagBREPCurve[...](int)
608de78e4feSLisandro Dalcin   ...
609de78e4feSLisandro Dalcin   // volumes
610de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
611de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
612de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
613de78e4feSLisandro Dalcin     numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int)
614de78e4feSLisandro Dalcin   ...
615de78e4feSLisandro Dalcin $EndEntities
616de78e4feSLisandro Dalcin */
6170598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh)
618de78e4feSLisandro Dalcin {
619698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
620698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
621698a79b9SLisandro Dalcin   long           index, num, lbuf[4];
622730356f6SLisandro Dalcin   int            dim, eid, numTags, *ibuf, t;
623698a79b9SLisandro Dalcin   PetscInt       count[4], i;
624a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
625de78e4feSLisandro Dalcin   PetscErrorCode ierr;
626de78e4feSLisandro Dalcin 
627de78e4feSLisandro Dalcin   PetscFunctionBegin;
628698a79b9SLisandro Dalcin   ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr);
629698a79b9SLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);}
630698a79b9SLisandro Dalcin   for (i = 0; i < 4; ++i) count[i] = lbuf[i];
6310598e1a2SLisandro Dalcin   ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr);
632de78e4feSLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
633730356f6SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
634730356f6SLisandro Dalcin       ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
635730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);}
6360598e1a2SLisandro Dalcin       ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr);
637de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
638de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);}
639de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
640de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
641698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr);
642de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
643730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
644de78e4feSLisandro Dalcin       entity->numTags = numTags = (int) PetscMin(num, 4);
645de78e4feSLisandro Dalcin       for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t];
646de78e4feSLisandro Dalcin       if (dim == 0) continue;
647de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
648de78e4feSLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);}
649698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr);
650de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr);
651730356f6SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);}
652de78e4feSLisandro Dalcin     }
653de78e4feSLisandro Dalcin   }
654de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
655de78e4feSLisandro Dalcin }
656de78e4feSLisandro Dalcin 
657de78e4feSLisandro Dalcin /*
658de78e4feSLisandro Dalcin $Nodes
659de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numNodes(unsigned long)
660de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long)
661de78e4feSLisandro Dalcin     tag(int) x(double) y(double) z(double)
662de78e4feSLisandro Dalcin     ...
663de78e4feSLisandro Dalcin   ...
664de78e4feSLisandro Dalcin $EndNodes
665de78e4feSLisandro Dalcin */
6660598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh)
667de78e4feSLisandro Dalcin {
668698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
669698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
6700598e1a2SLisandro Dalcin   long           block, node, n, numEntityBlocks, numTotalNodes, numNodes;
671de78e4feSLisandro Dalcin   int            info[3], nid;
6720598e1a2SLisandro Dalcin   GmshNodes      *nodes;
673de78e4feSLisandro Dalcin   PetscErrorCode ierr;
674de78e4feSLisandro Dalcin 
675de78e4feSLisandro Dalcin   PetscFunctionBegin;
676de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
677de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
678de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
679de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
6800598e1a2SLisandro Dalcin   ierr = GmshNodesCreate(numTotalNodes, &nodes);CHKERRQ(ierr);
6810598e1a2SLisandro Dalcin   mesh->numNodes = numTotalNodes;
6820598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
6830598e1a2SLisandro Dalcin   for (n = 0, block = 0; block < numEntityBlocks; ++block) {
684de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
685de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
686de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
687698a79b9SLisandro Dalcin     if (gmsh->binary) {
688277f51e8SBarry Smith       size_t nbytes = sizeof(int) + 3*sizeof(double);
689277f51e8SBarry Smith       char   *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */
690698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr);
691de78e4feSLisandro Dalcin       ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr);
6920598e1a2SLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++n) {
693de78e4feSLisandro Dalcin         char   *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int);
6940598e1a2SLisandro Dalcin         double *xyz = nodes->xyz + n*3;
69530815ce0SLisandro Dalcin         if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);}
69630815ce0SLisandro Dalcin         if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
697de78e4feSLisandro Dalcin         ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr);
698de78e4feSLisandro Dalcin         ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr);
699de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
700de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
7010598e1a2SLisandro Dalcin         nodes->id[n] = nid;
702de78e4feSLisandro Dalcin       }
703de78e4feSLisandro Dalcin     } else {
7040598e1a2SLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++n) {
7050598e1a2SLisandro Dalcin         double *xyz = nodes->xyz + n*3;
706de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
707de78e4feSLisandro Dalcin         ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
7080598e1a2SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);}
709de78e4feSLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);}
7100598e1a2SLisandro Dalcin         nodes->id[n] = nid;
711de78e4feSLisandro Dalcin       }
712de78e4feSLisandro Dalcin     }
713de78e4feSLisandro Dalcin   }
714de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
715de78e4feSLisandro Dalcin }
716de78e4feSLisandro Dalcin 
717de78e4feSLisandro Dalcin /*
718de78e4feSLisandro Dalcin $Elements
719de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numElements(unsigned long)
720de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long)
721de78e4feSLisandro Dalcin     tag(int) numVert[...](int)
722de78e4feSLisandro Dalcin     ...
723de78e4feSLisandro Dalcin   ...
724de78e4feSLisandro Dalcin $EndElements
725de78e4feSLisandro Dalcin */
7260598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh)
727de78e4feSLisandro Dalcin {
728698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
729698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
730de78e4feSLisandro Dalcin   long           c, block, numEntityBlocks, numTotalElements, elem, numElements;
731de78e4feSLisandro Dalcin   int            p, info[3], *ibuf = NULL;
7320598e1a2SLisandro Dalcin   int            eid, dim, cellType, numVerts, numNodes, numTags;
733a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
734de78e4feSLisandro Dalcin   GmshElement    *elements;
7350598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
736de78e4feSLisandro Dalcin   PetscErrorCode ierr;
737de78e4feSLisandro Dalcin 
738de78e4feSLisandro Dalcin   PetscFunctionBegin;
739de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
740de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);}
741de78e4feSLisandro Dalcin   ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
742de78e4feSLisandro Dalcin   if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);}
7430598e1a2SLisandro Dalcin   ierr = GmshElementsCreate(numTotalElements, &elements);CHKERRQ(ierr);
7440598e1a2SLisandro Dalcin   mesh->numElems = numTotalElements;
7450598e1a2SLisandro Dalcin   mesh->elements = elements;
746de78e4feSLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
747de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
748de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);}
749de78e4feSLisandro Dalcin     eid = info[0]; dim = info[1]; cellType = info[2];
7500598e1a2SLisandro Dalcin     ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr);
7510598e1a2SLisandro Dalcin     ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr);
7520598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
7530598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
754730356f6SLisandro Dalcin     numTags  = entity->numTags;
755de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
756de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);}
757698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr);
758de78e4feSLisandro Dalcin     ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr);
759de78e4feSLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);}
760de78e4feSLisandro Dalcin     for (elem = 0; elem < numElements; ++elem, ++c) {
761de78e4feSLisandro Dalcin       GmshElement *element = elements + c;
7620598e1a2SLisandro Dalcin       const int *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
7630598e1a2SLisandro Dalcin       element->id  = id[0];
764de78e4feSLisandro Dalcin       element->dim = dim;
765de78e4feSLisandro Dalcin       element->cellType = cellType;
7660598e1a2SLisandro Dalcin       element->numVerts = numVerts;
767de78e4feSLisandro Dalcin       element->numNodes = numNodes;
768de78e4feSLisandro Dalcin       element->numTags  = numTags;
7690598e1a2SLisandro Dalcin       ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr);
7700598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
7710598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = entity->tags[p];
772de78e4feSLisandro Dalcin     }
773de78e4feSLisandro Dalcin   }
774698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
775698a79b9SLisandro Dalcin }
776698a79b9SLisandro Dalcin 
7770598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[])
778698a79b9SLisandro Dalcin {
779698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
780698a79b9SLisandro Dalcin   int            fileFormat = gmsh->fileFormat;
781698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
782698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
783698a79b9SLisandro Dalcin   int            numPeriodic, snum, i;
784698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
7850598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
786698a79b9SLisandro Dalcin   PetscErrorCode ierr;
787698a79b9SLisandro Dalcin 
788698a79b9SLisandro Dalcin   PetscFunctionBegin;
789698a79b9SLisandro Dalcin   if (fileFormat == 22 || !binary) {
790698a79b9SLisandro Dalcin     ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
791698a79b9SLisandro Dalcin     snum = sscanf(line, "%d", &numPeriodic);
7920598e1a2SLisandro Dalcin     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
793698a79b9SLisandro Dalcin   } else {
794698a79b9SLisandro Dalcin     ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr);
795698a79b9SLisandro Dalcin     if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);}
796698a79b9SLisandro Dalcin   }
797698a79b9SLisandro Dalcin   for (i = 0; i < numPeriodic; i++) {
798698a79b9SLisandro Dalcin     int    ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode;
799698a79b9SLisandro Dalcin     long   j, nNodes;
800698a79b9SLisandro Dalcin     double affine[16];
801698a79b9SLisandro Dalcin 
802698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
803698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);
804698a79b9SLisandro Dalcin       snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag);
8050598e1a2SLisandro Dalcin       if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
806698a79b9SLisandro Dalcin     } else {
807698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr);
808698a79b9SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);}
809698a79b9SLisandro Dalcin       slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2];
810698a79b9SLisandro Dalcin     }
811698a79b9SLisandro Dalcin     (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */
812698a79b9SLisandro Dalcin 
813698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
814698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
815698a79b9SLisandro Dalcin       snum = sscanf(line, "%ld", &nNodes);
8164c500f23SPierre Jolivet       if (snum != 1) { /* discard transformation and try again */
817698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr);
818698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);
819698a79b9SLisandro Dalcin         snum = sscanf(line, "%ld", &nNodes);
8200598e1a2SLisandro Dalcin         if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
821698a79b9SLisandro Dalcin       }
822698a79b9SLisandro Dalcin     } else {
823698a79b9SLisandro Dalcin       ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
824698a79b9SLisandro Dalcin       if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
8254c500f23SPierre Jolivet       if (nNodes == -1) { /* discard transformation and try again */
826698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr);
827698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr);
828698a79b9SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);}
829698a79b9SLisandro Dalcin       }
830698a79b9SLisandro Dalcin     }
831698a79b9SLisandro Dalcin 
832698a79b9SLisandro Dalcin     for (j = 0; j < nNodes; j++) {
833698a79b9SLisandro Dalcin       if (fileFormat == 22 || !binary) {
834698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr);
835698a79b9SLisandro Dalcin         snum = sscanf(line, "%d %d", &slaveNode, &masterNode);
8360598e1a2SLisandro Dalcin         if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
837698a79b9SLisandro Dalcin       } else {
838698a79b9SLisandro Dalcin         ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr);
839698a79b9SLisandro Dalcin         if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);}
840698a79b9SLisandro Dalcin         slaveNode = ibuf[0]; masterNode = ibuf[1];
841698a79b9SLisandro Dalcin       }
8420598e1a2SLisandro Dalcin       slaveNode  = (int) nodeMap[slaveNode];
8430598e1a2SLisandro Dalcin       masterNode = (int) nodeMap[masterNode];
8440598e1a2SLisandro Dalcin       periodicMap[slaveNode] = masterNode;
845698a79b9SLisandro Dalcin     }
846698a79b9SLisandro Dalcin   }
847698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
848698a79b9SLisandro Dalcin }
849698a79b9SLisandro Dalcin 
8500598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
851698a79b9SLisandro Dalcin $Entities
852698a79b9SLisandro Dalcin   numPoints(size_t) numCurves(size_t)
853698a79b9SLisandro Dalcin     numSurfaces(size_t) numVolumes(size_t)
854698a79b9SLisandro Dalcin   pointTag(int) X(double) Y(double) Z(double)
855698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
856698a79b9SLisandro Dalcin   ...
857698a79b9SLisandro Dalcin   curveTag(int) minX(double) minY(double) minZ(double)
858698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
859698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
860698a79b9SLisandro Dalcin     numBoundingPoints(size_t) pointTag(int) ...
861698a79b9SLisandro Dalcin   ...
862698a79b9SLisandro Dalcin   surfaceTag(int) minX(double) minY(double) minZ(double)
863698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
864698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
865698a79b9SLisandro Dalcin     numBoundingCurves(size_t) curveTag(int) ...
866698a79b9SLisandro Dalcin   ...
867698a79b9SLisandro Dalcin   volumeTag(int) minX(double) minY(double) minZ(double)
868698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
869698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
870698a79b9SLisandro Dalcin     numBoundngSurfaces(size_t) surfaceTag(int) ...
871698a79b9SLisandro Dalcin   ...
872698a79b9SLisandro Dalcin $EndEntities
873698a79b9SLisandro Dalcin */
8740598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh)
875698a79b9SLisandro Dalcin {
876698a79b9SLisandro Dalcin   PetscInt       count[4], index, numTags, i;
877698a79b9SLisandro Dalcin   int            dim, eid, *tags = NULL;
878698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
879698a79b9SLisandro Dalcin   PetscErrorCode ierr;
880698a79b9SLisandro Dalcin 
881698a79b9SLisandro Dalcin   PetscFunctionBegin;
882698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr);
8830598e1a2SLisandro Dalcin   ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr);
884698a79b9SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
885698a79b9SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
886698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr);
8870598e1a2SLisandro Dalcin       ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr);
888698a79b9SLisandro Dalcin       ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr);
889698a79b9SLisandro Dalcin       ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr);
890698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr);
891698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr);
892698a79b9SLisandro Dalcin       entity->numTags = PetscMin(numTags, 4);
893698a79b9SLisandro Dalcin       for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i];
894698a79b9SLisandro Dalcin       if (dim == 0) continue;
895698a79b9SLisandro Dalcin       ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr);
896698a79b9SLisandro Dalcin       ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr);
897698a79b9SLisandro Dalcin       ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr);
898698a79b9SLisandro Dalcin     }
899698a79b9SLisandro Dalcin   }
900698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
901698a79b9SLisandro Dalcin }
902698a79b9SLisandro Dalcin 
90333a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
904698a79b9SLisandro Dalcin $Nodes
905698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numNodes(size_t)
906698a79b9SLisandro Dalcin     minNodeTag(size_t) maxNodeTag(size_t)
907698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t)
908698a79b9SLisandro Dalcin     nodeTag(size_t)
909698a79b9SLisandro Dalcin     ...
910698a79b9SLisandro Dalcin     x(double) y(double) z(double)
911698a79b9SLisandro Dalcin        < u(double; if parametric and entityDim = 1 or entityDim = 2) >
912698a79b9SLisandro Dalcin        < v(double; if parametric and entityDim = 2) >
913698a79b9SLisandro Dalcin     ...
914698a79b9SLisandro Dalcin   ...
915698a79b9SLisandro Dalcin $EndNodes
916698a79b9SLisandro Dalcin */
9170598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh)
918698a79b9SLisandro Dalcin {
919698a79b9SLisandro Dalcin   int            info[3];
9200598e1a2SLisandro Dalcin   PetscInt       sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, block, node;
9210598e1a2SLisandro Dalcin   GmshNodes      *nodes;
922698a79b9SLisandro Dalcin   PetscErrorCode ierr;
923698a79b9SLisandro Dalcin 
924698a79b9SLisandro Dalcin   PetscFunctionBegin;
925698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr);
9260598e1a2SLisandro Dalcin   numEntityBlocks = sizes[0]; numNodes = sizes[1];
9270598e1a2SLisandro Dalcin   ierr = GmshNodesCreate(numNodes, &nodes);CHKERRQ(ierr);
9280598e1a2SLisandro Dalcin   mesh->numNodes = numNodes;
9290598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
930698a79b9SLisandro Dalcin   for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) {
931698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
932698a79b9SLisandro Dalcin     if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported");
9330598e1a2SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr);
9340598e1a2SLisandro Dalcin     ierr = GmshReadSize(gmsh, nodes->id+node, numNodesBlock);CHKERRQ(ierr);
9350598e1a2SLisandro Dalcin     ierr = GmshReadDouble(gmsh, nodes->xyz+node*3, numNodesBlock*3);CHKERRQ(ierr);
936698a79b9SLisandro Dalcin   }
9370598e1a2SLisandro Dalcin   gmsh->nodeStart = sizes[2];
9380598e1a2SLisandro Dalcin   gmsh->nodeEnd   = sizes[3]+1;
939698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
940698a79b9SLisandro Dalcin }
941698a79b9SLisandro Dalcin 
94233a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
943698a79b9SLisandro Dalcin $Elements
944698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numElements(size_t)
945698a79b9SLisandro Dalcin     minElementTag(size_t) maxElementTag(size_t)
946698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t)
947698a79b9SLisandro Dalcin     elementTag(size_t) nodeTag(size_t) ...
948698a79b9SLisandro Dalcin     ...
949698a79b9SLisandro Dalcin   ...
950698a79b9SLisandro Dalcin $EndElements
951698a79b9SLisandro Dalcin */
9520598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh)
953698a79b9SLisandro Dalcin {
9540598e1a2SLisandro Dalcin   int            info[3], eid, dim, cellType;
9550598e1a2SLisandro Dalcin   PetscInt       sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p;
956698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
957698a79b9SLisandro Dalcin   GmshElement    *elements;
9580598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
959698a79b9SLisandro Dalcin   PetscErrorCode ierr;
960698a79b9SLisandro Dalcin 
961698a79b9SLisandro Dalcin   PetscFunctionBegin;
962698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr);
963698a79b9SLisandro Dalcin   numEntityBlocks = sizes[0]; numElements = sizes[1];
9640598e1a2SLisandro Dalcin   ierr = GmshElementsCreate(numElements, &elements);CHKERRQ(ierr);
9650598e1a2SLisandro Dalcin   mesh->numElems = numElements;
9660598e1a2SLisandro Dalcin   mesh->elements = elements;
967698a79b9SLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
968698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
969698a79b9SLisandro Dalcin     dim = info[0]; eid = info[1]; cellType = info[2];
9700598e1a2SLisandro Dalcin     ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr);
9710598e1a2SLisandro Dalcin     ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr);
9720598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
9730598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
974698a79b9SLisandro Dalcin     numTags  = entity->numTags;
975698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr);
976698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr);
977698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr);
978698a79b9SLisandro Dalcin     for (elem = 0; elem < numBlockElements; ++elem, ++c) {
979698a79b9SLisandro Dalcin       GmshElement *element = elements + c;
9800598e1a2SLisandro Dalcin       const PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
981698a79b9SLisandro Dalcin       element->id  = id[0];
982698a79b9SLisandro Dalcin       element->dim = dim;
983698a79b9SLisandro Dalcin       element->cellType = cellType;
9840598e1a2SLisandro Dalcin       element->numVerts = numVerts;
985698a79b9SLisandro Dalcin       element->numNodes = numNodes;
986698a79b9SLisandro Dalcin       element->numTags  = numTags;
9870598e1a2SLisandro Dalcin       ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr);
9880598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
9890598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = entity->tags[p];
990698a79b9SLisandro Dalcin     }
991698a79b9SLisandro Dalcin   }
992698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
993698a79b9SLisandro Dalcin }
994698a79b9SLisandro Dalcin 
9950598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
996698a79b9SLisandro Dalcin $Periodic
997698a79b9SLisandro Dalcin   numPeriodicLinks(size_t)
998698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) entityTagMaster(int)
999698a79b9SLisandro Dalcin   numAffine(size_t) value(double) ...
1000698a79b9SLisandro Dalcin   numCorrespondingNodes(size_t)
1001698a79b9SLisandro Dalcin     nodeTag(size_t) nodeTagMaster(size_t)
1002698a79b9SLisandro Dalcin     ...
1003698a79b9SLisandro Dalcin   ...
1004698a79b9SLisandro Dalcin $EndPeriodic
1005698a79b9SLisandro Dalcin */
10060598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[])
1007698a79b9SLisandro Dalcin {
1008698a79b9SLisandro Dalcin   int            info[3];
1009698a79b9SLisandro Dalcin   double         dbuf[16];
10100598e1a2SLisandro Dalcin   PetscInt       numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node;
10110598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
1012698a79b9SLisandro Dalcin   PetscErrorCode ierr;
1013698a79b9SLisandro Dalcin 
1014698a79b9SLisandro Dalcin   PetscFunctionBegin;
1015698a79b9SLisandro Dalcin   ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr);
1016698a79b9SLisandro Dalcin   for (link = 0; link < numPeriodicLinks; ++link) {
1017698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr);
1018698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr);
1019698a79b9SLisandro Dalcin     ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr);
1020698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr);
1021698a79b9SLisandro Dalcin     ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr);
1022698a79b9SLisandro Dalcin     ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr);
1023698a79b9SLisandro Dalcin     for (node = 0; node < numCorrespondingNodes; ++node) {
10240598e1a2SLisandro Dalcin       PetscInt slaveNode  = nodeMap[nodeTags[node*2+0]];
10250598e1a2SLisandro Dalcin       PetscInt masterNode = nodeMap[nodeTags[node*2+1]];
10260598e1a2SLisandro Dalcin       periodicMap[slaveNode] = masterNode;
1027698a79b9SLisandro Dalcin     }
1028698a79b9SLisandro Dalcin   }
1029698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
1030698a79b9SLisandro Dalcin }
1031698a79b9SLisandro Dalcin 
10320598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1033d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2
1034d6689ca9SLisandro Dalcin   version(ASCII double; currently 4.1)
1035d6689ca9SLisandro Dalcin   file-type(ASCII int; 0 for ASCII mode, 1 for binary mode)
1036d6689ca9SLisandro Dalcin   data-size(ASCII int; sizeof(size_t))
1037d6689ca9SLisandro Dalcin   < int with value one; only in binary mode, to detect endianness >
1038d6689ca9SLisandro Dalcin $EndMeshFormat
1039d6689ca9SLisandro Dalcin */
10400598e1a2SLisandro Dalcin static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh)
1041698a79b9SLisandro Dalcin {
1042698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
1043698a79b9SLisandro Dalcin   int            snum, fileType, fileFormat, dataSize, checkEndian;
1044698a79b9SLisandro Dalcin   float          version;
1045698a79b9SLisandro Dalcin   PetscErrorCode ierr;
1046698a79b9SLisandro Dalcin 
1047698a79b9SLisandro Dalcin   PetscFunctionBegin;
1048698a79b9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr);
1049698a79b9SLisandro Dalcin   snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
10500598e1a2SLisandro Dalcin   if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
10510598e1a2SLisandro Dalcin   if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
10520598e1a2SLisandro Dalcin   if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
10530598e1a2SLisandro Dalcin   if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
10540598e1a2SLisandro Dalcin   if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII");
10550598e1a2SLisandro Dalcin   if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary");
1056698a79b9SLisandro Dalcin   fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */
10570598e1a2SLisandro Dalcin   if (fileFormat <= 40 && dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize);
10580598e1a2SLisandro Dalcin   if (fileFormat >= 41 && dataSize != sizeof(int) && dataSize != sizeof(PetscInt64)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize);
1059698a79b9SLisandro Dalcin   gmsh->fileFormat = fileFormat;
1060698a79b9SLisandro Dalcin   gmsh->dataSize = dataSize;
1061698a79b9SLisandro Dalcin   gmsh->byteSwap = PETSC_FALSE;
1062698a79b9SLisandro Dalcin   if (gmsh->binary) {
1063698a79b9SLisandro Dalcin     ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr);
1064698a79b9SLisandro Dalcin     if (checkEndian != 1) {
1065698a79b9SLisandro Dalcin       ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr);
10660598e1a2SLisandro Dalcin       if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line);
1067698a79b9SLisandro Dalcin       gmsh->byteSwap = PETSC_TRUE;
1068698a79b9SLisandro Dalcin     }
1069698a79b9SLisandro Dalcin   }
1070698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
1071698a79b9SLisandro Dalcin }
1072698a79b9SLisandro Dalcin 
10731f643af3SLisandro Dalcin /*
10741f643af3SLisandro Dalcin PhysicalNames
10751f643af3SLisandro Dalcin   numPhysicalNames(ASCII int)
10761f643af3SLisandro Dalcin   dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max)
10771f643af3SLisandro Dalcin   ...
10781f643af3SLisandro Dalcin $EndPhysicalNames
10791f643af3SLisandro Dalcin */
10800598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh)
1081698a79b9SLisandro Dalcin {
10821f643af3SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q;
10831f643af3SLisandro Dalcin   int            snum, numRegions, region, dim, tag;
1084698a79b9SLisandro Dalcin   PetscErrorCode ierr;
1085698a79b9SLisandro Dalcin 
1086698a79b9SLisandro Dalcin   PetscFunctionBegin;
1087698a79b9SLisandro Dalcin   ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr);
1088698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &numRegions);
10890598e1a2SLisandro Dalcin   if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1090698a79b9SLisandro Dalcin   for (region = 0; region < numRegions; ++region) {
10911f643af3SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr);
10921f643af3SLisandro Dalcin     snum = sscanf(line, "%d %d", &dim, &tag);
10930598e1a2SLisandro Dalcin     if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
10941f643af3SLisandro Dalcin     ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr);
10951f643af3SLisandro Dalcin     ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr);
10960598e1a2SLisandro Dalcin     if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
10971f643af3SLisandro Dalcin     ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr);
10980598e1a2SLisandro Dalcin     if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
10991f643af3SLisandro Dalcin     ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr);
1100698a79b9SLisandro Dalcin   }
1101de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
1102de78e4feSLisandro Dalcin }
1103de78e4feSLisandro Dalcin 
11040598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh)
110596ca5757SLisandro Dalcin {
11060598e1a2SLisandro Dalcin   PetscErrorCode ierr;
11070598e1a2SLisandro Dalcin 
11080598e1a2SLisandro Dalcin   PetscFunctionBegin;
11090598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11100598e1a2SLisandro Dalcin   case 41: ierr = GmshReadEntities_v41(gmsh, mesh);CHKERRQ(ierr); break;
11110598e1a2SLisandro Dalcin   default: ierr = GmshReadEntities_v40(gmsh, mesh);CHKERRQ(ierr); break;
111296ca5757SLisandro Dalcin   }
11130598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
11140598e1a2SLisandro Dalcin }
11150598e1a2SLisandro Dalcin 
11160598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh)
11170598e1a2SLisandro Dalcin {
11180598e1a2SLisandro Dalcin   PetscErrorCode ierr;
11190598e1a2SLisandro Dalcin 
11200598e1a2SLisandro Dalcin   PetscFunctionBegin;
11210598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11220598e1a2SLisandro Dalcin   case 41: ierr = GmshReadNodes_v41(gmsh, mesh);CHKERRQ(ierr); break;
11230598e1a2SLisandro Dalcin   case 40: ierr = GmshReadNodes_v40(gmsh, mesh);CHKERRQ(ierr); break;
11240598e1a2SLisandro Dalcin   default: ierr = GmshReadNodes_v22(gmsh, mesh);CHKERRQ(ierr); break;
11250598e1a2SLisandro Dalcin   }
11260598e1a2SLisandro Dalcin 
11270598e1a2SLisandro Dalcin   { /* Gmsh v2.2/v4.0 does not provide min/max node tags */
11280598e1a2SLisandro Dalcin     if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) {
11290598e1a2SLisandro Dalcin       PetscInt  tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n;
11300598e1a2SLisandro Dalcin       GmshNodes *nodes = mesh->nodelist;
11310598e1a2SLisandro Dalcin       for (n = 0; n < mesh->numNodes; ++n) {
11320598e1a2SLisandro Dalcin         const PetscInt tag = nodes->id[n];
11330598e1a2SLisandro Dalcin         tagMin = PetscMin(tag, tagMin);
11340598e1a2SLisandro Dalcin         tagMax = PetscMax(tag, tagMax);
11350598e1a2SLisandro Dalcin       }
11360598e1a2SLisandro Dalcin       gmsh->nodeStart = tagMin;
11370598e1a2SLisandro Dalcin       gmsh->nodeEnd   = tagMax+1;
11380598e1a2SLisandro Dalcin     }
11390598e1a2SLisandro Dalcin   }
11400598e1a2SLisandro Dalcin 
11410598e1a2SLisandro Dalcin   { /* Support for sparse node tags */
11420598e1a2SLisandro Dalcin     PetscInt  n, t;
11430598e1a2SLisandro Dalcin     GmshNodes *nodes = mesh->nodelist;
11440598e1a2SLisandro Dalcin     ierr = PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf);CHKERRQ(ierr);
11450598e1a2SLisandro Dalcin     for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT;
11460598e1a2SLisandro Dalcin     gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart;
11470598e1a2SLisandro Dalcin     for (n = 0; n < mesh->numNodes; ++n) {
11480598e1a2SLisandro Dalcin       const PetscInt tag = nodes->id[n];
11490598e1a2SLisandro Dalcin       if (gmsh->nodeMap[tag] >= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %D", tag);
11500598e1a2SLisandro Dalcin       gmsh->nodeMap[tag] = n;
11510598e1a2SLisandro Dalcin     }
11520598e1a2SLisandro Dalcin   }
11530598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
11540598e1a2SLisandro Dalcin }
11550598e1a2SLisandro Dalcin 
11560598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh)
11570598e1a2SLisandro Dalcin {
11580598e1a2SLisandro Dalcin   PetscErrorCode ierr;
11590598e1a2SLisandro Dalcin 
11600598e1a2SLisandro Dalcin   PetscFunctionBegin;
11610598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11620598e1a2SLisandro Dalcin   case 41: ierr = GmshReadElements_v41(gmsh, mesh);CHKERRQ(ierr); break;
11630598e1a2SLisandro Dalcin   case 40: ierr = GmshReadElements_v40(gmsh, mesh);CHKERRQ(ierr); break;
11640598e1a2SLisandro Dalcin   default: ierr = GmshReadElements_v22(gmsh, mesh);CHKERRQ(ierr); break;
11650598e1a2SLisandro Dalcin   }
11660598e1a2SLisandro Dalcin 
11670598e1a2SLisandro Dalcin   { /* Reorder elements by codimension and polytope type */
11680598e1a2SLisandro Dalcin     PetscInt    ne = mesh->numElems;
11690598e1a2SLisandro Dalcin     GmshElement *elements = mesh->elements;
1170*066ea43fSLisandro Dalcin     PetscInt    keymap[GMSH_NUM_POLYTOPES], nk = 0;
1171*066ea43fSLisandro Dalcin     PetscInt    offset[GMSH_NUM_POLYTOPES+1], e, k;
11720598e1a2SLisandro Dalcin 
1173*066ea43fSLisandro Dalcin     for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT;
11740598e1a2SLisandro Dalcin     ierr = PetscMemzero(offset,sizeof(offset));CHKERRQ(ierr);
11750598e1a2SLisandro Dalcin 
1176*066ea43fSLisandro Dalcin     keymap[GMSH_TET] = nk++;
1177*066ea43fSLisandro Dalcin     keymap[GMSH_HEX] = nk++;
1178*066ea43fSLisandro Dalcin     keymap[GMSH_PRI] = nk++;
1179*066ea43fSLisandro Dalcin     keymap[GMSH_PYR] = nk++;
1180*066ea43fSLisandro Dalcin     keymap[GMSH_TRI] = nk++;
1181*066ea43fSLisandro Dalcin     keymap[GMSH_QUA] = nk++;
1182*066ea43fSLisandro Dalcin     keymap[GMSH_SEG] = nk++;
1183*066ea43fSLisandro Dalcin     keymap[GMSH_VTX] = nk++;
11840598e1a2SLisandro Dalcin 
11850598e1a2SLisandro Dalcin     ierr = GmshElementsCreate(mesh->numElems, &mesh->elements);CHKERRQ(ierr);
11860598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope]
11870598e1a2SLisandro Dalcin     for (e = 0; e < ne; ++e) offset[1+key(e)]++;
11880598e1a2SLisandro Dalcin     for (k = 1; k < nk; ++k) offset[k] += offset[k-1];
11890598e1a2SLisandro Dalcin     for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e];
11900598e1a2SLisandro Dalcin #undef key
11910598e1a2SLisandro Dalcin     ierr = GmshElementsDestroy(&elements);CHKERRQ(ierr);
1192*066ea43fSLisandro Dalcin   }
11930598e1a2SLisandro Dalcin 
1194*066ea43fSLisandro Dalcin   { /* Mesh dimension and order */
1195*066ea43fSLisandro Dalcin     GmshElement *elem = mesh->numElems ? mesh->elements : NULL;
1196*066ea43fSLisandro Dalcin     mesh->dim   = elem ? GmshCellMap[elem->cellType].dim   : 0;
1197*066ea43fSLisandro Dalcin     mesh->order = elem ? GmshCellMap[elem->cellType].order : 0;
11980598e1a2SLisandro Dalcin   }
11990598e1a2SLisandro Dalcin 
12000598e1a2SLisandro Dalcin   {
12010598e1a2SLisandro Dalcin     PetscBT  vtx;
12020598e1a2SLisandro Dalcin     PetscInt dim = mesh->dim, e, n, v;
12030598e1a2SLisandro Dalcin 
12040598e1a2SLisandro Dalcin     ierr = PetscBTCreate(mesh->numNodes, &vtx);CHKERRQ(ierr);
12050598e1a2SLisandro Dalcin 
12060598e1a2SLisandro Dalcin     /* Compute number of cells and set of vertices */
12070598e1a2SLisandro Dalcin     mesh->numCells = 0;
12080598e1a2SLisandro Dalcin     for (e = 0; e < mesh->numElems; ++e) {
12090598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + e;
12100598e1a2SLisandro Dalcin       if (elem->dim == dim && dim > 0) mesh->numCells++;
12110598e1a2SLisandro Dalcin       for (v = 0; v < elem->numVerts; v++) {
12120598e1a2SLisandro Dalcin         ierr = PetscBTSet(vtx, elem->nodes[v]);CHKERRQ(ierr);
12130598e1a2SLisandro Dalcin       }
12140598e1a2SLisandro Dalcin     }
12150598e1a2SLisandro Dalcin 
12160598e1a2SLisandro Dalcin     /* Compute numbering for vertices */
12170598e1a2SLisandro Dalcin     mesh->numVerts = 0;
12180598e1a2SLisandro Dalcin     ierr = PetscMalloc1(mesh->numNodes, &mesh->vertexMap);CHKERRQ(ierr);
12190598e1a2SLisandro Dalcin     for (n = 0; n < mesh->numNodes; ++n)
12200598e1a2SLisandro Dalcin       mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT;
12210598e1a2SLisandro Dalcin 
12220598e1a2SLisandro Dalcin     ierr = PetscBTDestroy(&vtx);CHKERRQ(ierr);
12230598e1a2SLisandro Dalcin   }
12240598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
12250598e1a2SLisandro Dalcin }
12260598e1a2SLisandro Dalcin 
12270598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh)
12280598e1a2SLisandro Dalcin {
12290598e1a2SLisandro Dalcin   PetscInt       n;
12300598e1a2SLisandro Dalcin   PetscErrorCode ierr;
12310598e1a2SLisandro Dalcin 
12320598e1a2SLisandro Dalcin   PetscFunctionBegin;
12330598e1a2SLisandro Dalcin   ierr = PetscMalloc1(mesh->numNodes, &mesh->periodMap);CHKERRQ(ierr);
12340598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n;
12350598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
12360598e1a2SLisandro Dalcin   case 41: ierr = GmshReadPeriodic_v41(gmsh, mesh->periodMap);CHKERRQ(ierr); break;
12370598e1a2SLisandro Dalcin   default: ierr = GmshReadPeriodic_v40(gmsh, mesh->periodMap);CHKERRQ(ierr); break;
12380598e1a2SLisandro Dalcin   }
12390598e1a2SLisandro Dalcin 
12400598e1a2SLisandro Dalcin   /* Find canonical master nodes */
12410598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12420598e1a2SLisandro Dalcin     while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]])
12430598e1a2SLisandro Dalcin       mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]];
12440598e1a2SLisandro Dalcin 
12450598e1a2SLisandro Dalcin   /* Renumber vertices (filter out slaves) */
12460598e1a2SLisandro Dalcin   mesh->numVerts = 0;
12470598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12480598e1a2SLisandro Dalcin     if (mesh->vertexMap[n] >= 0)   /* is vertex */
12490598e1a2SLisandro Dalcin       if (mesh->periodMap[n] == n) /* is master */
12500598e1a2SLisandro Dalcin         mesh->vertexMap[n] = mesh->numVerts++;
12510598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12520598e1a2SLisandro Dalcin     if (mesh->vertexMap[n] >= 0)   /* is vertex */
12530598e1a2SLisandro Dalcin       if (mesh->periodMap[n] != n) /* is slave  */
12540598e1a2SLisandro Dalcin         mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]];
12550598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
12560598e1a2SLisandro Dalcin }
12570598e1a2SLisandro Dalcin 
1258*066ea43fSLisandro Dalcin #define DM_POLYTOPE_VERTEX  DM_POLYTOPE_POINT
1259*066ea43fSLisandro Dalcin #define DM_POLYTOPE_PYRAMID DM_POLYTOPE_UNKNOWN
1260*066ea43fSLisandro Dalcin static const DMPolytopeType DMPolytopeMap[] = {
1261*066ea43fSLisandro Dalcin   /* GMSH_VTX */ DM_POLYTOPE_VERTEX,
1262*066ea43fSLisandro Dalcin   /* GMSH_SEG */ DM_POLYTOPE_SEGMENT,
1263*066ea43fSLisandro Dalcin   /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE,
1264*066ea43fSLisandro Dalcin   /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL,
1265*066ea43fSLisandro Dalcin   /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON,
1266*066ea43fSLisandro Dalcin   /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON,
1267*066ea43fSLisandro Dalcin   /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM,
1268*066ea43fSLisandro Dalcin   /* GMSH_PYR */ DM_POLYTOPE_PYRAMID,
1269*066ea43fSLisandro Dalcin   DM_POLYTOPE_UNKNOWN
1270*066ea43fSLisandro Dalcin };
12710598e1a2SLisandro Dalcin 
12720598e1a2SLisandro Dalcin PETSC_STATIC_INLINE DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType)
12730598e1a2SLisandro Dalcin {
1274*066ea43fSLisandro Dalcin   return DMPolytopeMap[GmshCellMap[cellType].polytope];
1275*066ea43fSLisandro Dalcin }
1276*066ea43fSLisandro Dalcin 
1277*066ea43fSLisandro Dalcin static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem)
1278*066ea43fSLisandro Dalcin {
1279*066ea43fSLisandro Dalcin   DM              K;
1280*066ea43fSLisandro Dalcin   PetscSpace      P;
1281*066ea43fSLisandro Dalcin   PetscDualSpace  Q;
1282*066ea43fSLisandro Dalcin   PetscQuadrature q, fq;
1283*066ea43fSLisandro Dalcin   PetscBool       isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1284*066ea43fSLisandro Dalcin   PetscBool       endpoint = PETSC_TRUE;
1285*066ea43fSLisandro Dalcin   char            name[32];
1286*066ea43fSLisandro Dalcin   PetscErrorCode  ierr;
1287*066ea43fSLisandro Dalcin 
1288*066ea43fSLisandro Dalcin   PetscFunctionBegin;
1289*066ea43fSLisandro Dalcin   /* Create space */
1290*066ea43fSLisandro Dalcin   ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr);
1291*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
1292*066ea43fSLisandro Dalcin   ierr = PetscSpacePolynomialSetTensor(P, isTensor);CHKERRQ(ierr);
1293*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr);
1294*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr);
1295*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetDegree(P, k, PETSC_DETERMINE);CHKERRQ(ierr);
1296*066ea43fSLisandro Dalcin   if (prefix) {
1297*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr);
1298*066ea43fSLisandro Dalcin     ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr);
1299*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) P, NULL);CHKERRQ(ierr);
1300*066ea43fSLisandro Dalcin     ierr = PetscSpaceGetDegree(P, &k, NULL);CHKERRQ(ierr);
1301*066ea43fSLisandro Dalcin   }
1302*066ea43fSLisandro Dalcin   ierr = PetscSpaceSetUp(P);CHKERRQ(ierr);
1303*066ea43fSLisandro Dalcin   /* Create dual space */
1304*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr);
1305*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr);
1306*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetTensor(Q, isTensor);CHKERRQ(ierr);
1307*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetContinuity(Q, continuity);CHKERRQ(ierr);
1308*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0);CHKERRQ(ierr);
1309*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr);
1310*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetOrder(Q, k);CHKERRQ(ierr);
1311*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr);
1312*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr);
1313*066ea43fSLisandro Dalcin   ierr = DMDestroy(&K);CHKERRQ(ierr);
1314*066ea43fSLisandro Dalcin   if (prefix) {
1315*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr);
1316*066ea43fSLisandro Dalcin     ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr);
1317*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, NULL);CHKERRQ(ierr);
1318*066ea43fSLisandro Dalcin   }
1319*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr);
1320*066ea43fSLisandro Dalcin   /* Create quadrature */
1321*066ea43fSLisandro Dalcin   if (isSimplex) {
1322*066ea43fSLisandro Dalcin     ierr = PetscDTStroudConicalQuadrature(dim,   1, k+1, -1, +1, &q);CHKERRQ(ierr);
1323*066ea43fSLisandro Dalcin     ierr = PetscDTStroudConicalQuadrature(dim-1, 1, k+1, -1, +1, &fq);CHKERRQ(ierr);
1324*066ea43fSLisandro Dalcin   } else {
1325*066ea43fSLisandro Dalcin     ierr = PetscDTGaussTensorQuadrature(dim,   1, k+1, -1, +1, &q);CHKERRQ(ierr);
1326*066ea43fSLisandro Dalcin     ierr = PetscDTGaussTensorQuadrature(dim-1, 1, k+1, -1, +1, &fq);CHKERRQ(ierr);
1327*066ea43fSLisandro Dalcin   }
1328*066ea43fSLisandro Dalcin   /* Create finite element */
1329*066ea43fSLisandro Dalcin   ierr = PetscFECreate(comm, fem);CHKERRQ(ierr);
1330*066ea43fSLisandro Dalcin   ierr = PetscFESetType(*fem, PETSCFEBASIC);CHKERRQ(ierr);
1331*066ea43fSLisandro Dalcin   ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr);
1332*066ea43fSLisandro Dalcin   ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr);
1333*066ea43fSLisandro Dalcin   ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr);
1334*066ea43fSLisandro Dalcin   ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr);
1335*066ea43fSLisandro Dalcin   ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr);
1336*066ea43fSLisandro Dalcin   if (prefix) {
1337*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr);
1338*066ea43fSLisandro Dalcin     ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr);
1339*066ea43fSLisandro Dalcin     ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, NULL);CHKERRQ(ierr);
1340*066ea43fSLisandro Dalcin   }
1341*066ea43fSLisandro Dalcin   ierr = PetscFESetUp(*fem);CHKERRQ(ierr);
1342*066ea43fSLisandro Dalcin   /* Cleanup */
1343*066ea43fSLisandro Dalcin   ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr);
1344*066ea43fSLisandro Dalcin   ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr);
1345*066ea43fSLisandro Dalcin   ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr);
1346*066ea43fSLisandro Dalcin   ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr);
1347*066ea43fSLisandro Dalcin   /* Set finite element name */
1348*066ea43fSLisandro Dalcin   ierr = PetscSNPrintf(name, sizeof(name), "%s%D", isSimplex? "P" : "Q", k);CHKERRQ(ierr);
1349*066ea43fSLisandro Dalcin   ierr = PetscFESetName(*fem, name);CHKERRQ(ierr);
1350*066ea43fSLisandro Dalcin   PetscFunctionReturn(0);
135196ca5757SLisandro Dalcin }
135296ca5757SLisandro Dalcin 
1353d6689ca9SLisandro Dalcin /*@C
1354d6689ca9SLisandro Dalcin   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
1355d6689ca9SLisandro Dalcin 
1356d6689ca9SLisandro Dalcin + comm        - The MPI communicator
1357d6689ca9SLisandro Dalcin . filename    - Name of the Gmsh file
1358d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh
1359d6689ca9SLisandro Dalcin 
1360d6689ca9SLisandro Dalcin   Output Parameter:
1361d6689ca9SLisandro Dalcin . dm  - The DM object representing the mesh
1362d6689ca9SLisandro Dalcin 
1363d6689ca9SLisandro Dalcin   Level: beginner
1364d6689ca9SLisandro Dalcin 
1365d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
1366d6689ca9SLisandro Dalcin @*/
1367d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1368d6689ca9SLisandro Dalcin {
1369d6689ca9SLisandro Dalcin   PetscViewer     viewer;
1370d6689ca9SLisandro Dalcin   PetscMPIInt     rank;
1371d6689ca9SLisandro Dalcin   int             fileType;
1372d6689ca9SLisandro Dalcin   PetscViewerType vtype;
1373d6689ca9SLisandro Dalcin   PetscErrorCode  ierr;
1374d6689ca9SLisandro Dalcin 
1375d6689ca9SLisandro Dalcin   PetscFunctionBegin;
1376d6689ca9SLisandro Dalcin   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1377d6689ca9SLisandro Dalcin 
1378d6689ca9SLisandro Dalcin   /* Determine Gmsh file type (ASCII or binary) from file header */
1379d6689ca9SLisandro Dalcin   if (!rank) {
13800598e1a2SLisandro Dalcin     GmshFile    gmsh[1];
1381d6689ca9SLisandro Dalcin     char        line[PETSC_MAX_PATH_LEN];
1382d6689ca9SLisandro Dalcin     int         snum;
1383d6689ca9SLisandro Dalcin     float       version;
1384d6689ca9SLisandro Dalcin 
1385580bdb30SBarry Smith     ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr);
1386d6689ca9SLisandro Dalcin     ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr);
1387d6689ca9SLisandro Dalcin     ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
1388d6689ca9SLisandro Dalcin     ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr);
1389d6689ca9SLisandro Dalcin     ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr);
1390d6689ca9SLisandro Dalcin     /* Read only the first two lines of the Gmsh file */
1391d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1392d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr);
1393d6689ca9SLisandro Dalcin     ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr);
1394d6689ca9SLisandro Dalcin     snum = sscanf(line, "%f %d", &version, &fileType);
13950598e1a2SLisandro Dalcin     if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
13960598e1a2SLisandro Dalcin     if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
13970598e1a2SLisandro Dalcin     if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
13980598e1a2SLisandro Dalcin     if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
1399d6689ca9SLisandro Dalcin     ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr);
1400d6689ca9SLisandro Dalcin   }
1401d6689ca9SLisandro Dalcin   ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr);
1402d6689ca9SLisandro Dalcin   vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
1403d6689ca9SLisandro Dalcin 
1404d6689ca9SLisandro Dalcin   /* Create appropriate viewer and build plex */
1405d6689ca9SLisandro Dalcin   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
1406d6689ca9SLisandro Dalcin   ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr);
1407d6689ca9SLisandro Dalcin   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
1408d6689ca9SLisandro Dalcin   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
1409d6689ca9SLisandro Dalcin   ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr);
1410d6689ca9SLisandro Dalcin   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1411d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
1412d6689ca9SLisandro Dalcin }
1413d6689ca9SLisandro Dalcin 
1414331830f3SMatthew G. Knepley /*@
14157d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
1416331830f3SMatthew G. Knepley 
1417d083f849SBarry Smith   Collective
1418331830f3SMatthew G. Knepley 
1419331830f3SMatthew G. Knepley   Input Parameters:
1420331830f3SMatthew G. Knepley + comm  - The MPI communicator
1421331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
1422331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
1423331830f3SMatthew G. Knepley 
1424331830f3SMatthew G. Knepley   Output Parameter:
1425331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
1426331830f3SMatthew G. Knepley 
1427698a79b9SLisandro Dalcin   Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format
1428331830f3SMatthew G. Knepley 
1429331830f3SMatthew G. Knepley   Level: beginner
1430331830f3SMatthew G. Knepley 
1431331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
1432331830f3SMatthew G. Knepley @*/
1433331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
1434331830f3SMatthew G. Knepley {
14350598e1a2SLisandro Dalcin   GmshMesh      *mesh = NULL;
143611c56cb0SLisandro Dalcin   PetscViewer    parentviewer = NULL;
14370598e1a2SLisandro Dalcin   PetscBT        periodicVerts = NULL;
14380598e1a2SLisandro Dalcin   PetscBT        periodicCells = NULL;
1439*066ea43fSLisandro Dalcin   DM             cdm;
1440331830f3SMatthew G. Knepley   PetscSection   coordSection;
1441331830f3SMatthew G. Knepley   Vec            coordinates;
1442*066ea43fSLisandro Dalcin   PetscInt       dim = 0, coordDim = -1, order = 0;
14430598e1a2SLisandro Dalcin   PetscInt       numNodes = 0, numElems = 0, numVerts = 0, numCells = 0;
1444*066ea43fSLisandro Dalcin   PetscInt       cell, cone[8], e, n, v, d;
14450598e1a2SLisandro Dalcin   PetscBool      binary, usemarker = PETSC_FALSE;
14460598e1a2SLisandro Dalcin   PetscBool      hybrid = interpolate, periodic = PETSC_TRUE;
1447*066ea43fSLisandro Dalcin   PetscBool      highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE;
1448*066ea43fSLisandro Dalcin   PetscBool      isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE;
14490598e1a2SLisandro Dalcin   PetscMPIInt    rank;
1450331830f3SMatthew G. Knepley   PetscErrorCode ierr;
1451331830f3SMatthew G. Knepley 
1452331830f3SMatthew G. Knepley   PetscFunctionBegin;
1453331830f3SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1454ef12879bSLisandro Dalcin   ierr = PetscObjectOptionsBegin((PetscObject)viewer);CHKERRQ(ierr);
1455ef12879bSLisandro Dalcin   ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Gmsh options");CHKERRQ(ierr);
14560598e1a2SLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", hybrid, &hybrid, NULL);CHKERRQ(ierr);
1457ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL);CHKERRQ(ierr);
1458*066ea43fSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_highorder","Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet);CHKERRQ(ierr);
1459*066ea43fSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL);CHKERRQ(ierr);
1460ef12879bSLisandro Dalcin   ierr = PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL);CHKERRQ(ierr);
14610598e1a2SLisandro Dalcin   ierr = PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE);CHKERRQ(ierr);
1462ef12879bSLisandro Dalcin   ierr = PetscOptionsTail();CHKERRQ(ierr);
1463ef12879bSLisandro Dalcin   ierr = PetscOptionsEnd();CHKERRQ(ierr);
14640598e1a2SLisandro Dalcin 
14650598e1a2SLisandro Dalcin   ierr = GmshCellInfoSetUp();CHKERRQ(ierr);
146611c56cb0SLisandro Dalcin 
1467331830f3SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1468331830f3SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
14690598e1a2SLisandro Dalcin   ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr);
147011c56cb0SLisandro Dalcin 
147111c56cb0SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr);
147211c56cb0SLisandro Dalcin 
147311c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
14743b46f529SLisandro Dalcin   if (binary) {
147511c56cb0SLisandro Dalcin     parentviewer = viewer;
147611c56cb0SLisandro Dalcin     ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
147711c56cb0SLisandro Dalcin   }
147811c56cb0SLisandro Dalcin 
14793b46f529SLisandro Dalcin   if (!rank) {
14800598e1a2SLisandro Dalcin     GmshFile  gmsh[1];
1481698a79b9SLisandro Dalcin     char      line[PETSC_MAX_PATH_LEN];
1482698a79b9SLisandro Dalcin     PetscBool match;
1483331830f3SMatthew G. Knepley 
1484580bdb30SBarry Smith     ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr);
1485698a79b9SLisandro Dalcin     gmsh->viewer = viewer;
1486698a79b9SLisandro Dalcin     gmsh->binary = binary;
1487698a79b9SLisandro Dalcin 
14880598e1a2SLisandro Dalcin     ierr = GmshMeshCreate(&mesh);CHKERRQ(ierr);
14890598e1a2SLisandro Dalcin 
1490698a79b9SLisandro Dalcin     /* Read mesh format */
1491d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1492d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr);
14930598e1a2SLisandro Dalcin     ierr = GmshReadMeshFormat(gmsh);CHKERRQ(ierr);
1494d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr);
149511c56cb0SLisandro Dalcin 
1496dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
1497d6689ca9SLisandro Dalcin     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1498d6689ca9SLisandro Dalcin     ierr = GmshMatch(gmsh, "$PhysicalNames", line, &match);CHKERRQ(ierr);
1499dc0ede3bSMatthew G. Knepley     if (match) {
15000598e1a2SLisandro Dalcin       ierr = GmshExpect(gmsh, "$PhysicalNames", line);CHKERRQ(ierr);
15010598e1a2SLisandro Dalcin       ierr = GmshReadPhysicalNames(gmsh);CHKERRQ(ierr);
1502d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr);
1503698a79b9SLisandro Dalcin       /* Initial read for entity section */
1504d6689ca9SLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1505dc0ede3bSMatthew G. Knepley     }
150611c56cb0SLisandro Dalcin 
1507de78e4feSLisandro Dalcin     /* Read entities */
1508698a79b9SLisandro Dalcin     if (gmsh->fileFormat >= 40) {
1509d6689ca9SLisandro Dalcin       ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr);
15100598e1a2SLisandro Dalcin       ierr = GmshReadEntities(gmsh, mesh);CHKERRQ(ierr);
1511d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr);
1512698a79b9SLisandro Dalcin       /* Initial read for nodes section */
1513d6689ca9SLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1514de78e4feSLisandro Dalcin     }
1515de78e4feSLisandro Dalcin 
1516698a79b9SLisandro Dalcin     /* Read nodes */
1517d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr);
15180598e1a2SLisandro Dalcin     ierr = GmshReadNodes(gmsh, mesh);CHKERRQ(ierr);
1519d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr);
152011c56cb0SLisandro Dalcin 
1521698a79b9SLisandro Dalcin     /* Read elements */
1522feb237baSPierre Jolivet     ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1523d6689ca9SLisandro Dalcin     ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr);
15240598e1a2SLisandro Dalcin     ierr = GmshReadElements(gmsh, mesh);CHKERRQ(ierr);
1525d6689ca9SLisandro Dalcin     ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr);
1526de78e4feSLisandro Dalcin 
15270598e1a2SLisandro Dalcin     /* Read periodic section (OPTIONAL) */
1528698a79b9SLisandro Dalcin     if (periodic) {
1529ef12879bSLisandro Dalcin       ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr);
1530ef12879bSLisandro Dalcin       ierr = GmshMatch(gmsh, "$Periodic", line, &periodic);CHKERRQ(ierr);
1531ef12879bSLisandro Dalcin     }
1532ef12879bSLisandro Dalcin     if (periodic) {
1533d6689ca9SLisandro Dalcin       ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr);
15340598e1a2SLisandro Dalcin       ierr = GmshReadPeriodic(gmsh, mesh);CHKERRQ(ierr);
1535d6689ca9SLisandro Dalcin       ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr);
1536698a79b9SLisandro Dalcin     }
1537698a79b9SLisandro Dalcin 
1538698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr);
1539698a79b9SLisandro Dalcin     ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr);
15400598e1a2SLisandro Dalcin     ierr = PetscFree(gmsh->nbuf);CHKERRQ(ierr);
15410598e1a2SLisandro Dalcin 
15420598e1a2SLisandro Dalcin     dim       = mesh->dim;
1543*066ea43fSLisandro Dalcin     order     = mesh->order;
15440598e1a2SLisandro Dalcin     numNodes  = mesh->numNodes;
15450598e1a2SLisandro Dalcin     numElems  = mesh->numElems;
15460598e1a2SLisandro Dalcin     numVerts  = mesh->numVerts;
15470598e1a2SLisandro Dalcin     numCells  = mesh->numCells;
1548*066ea43fSLisandro Dalcin 
1549*066ea43fSLisandro Dalcin     {
1550*066ea43fSLisandro Dalcin       GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL;
1551*066ea43fSLisandro Dalcin       GmshElement *elemB = elemA ? elemA + mesh->numCells - 1  : NULL;
1552*066ea43fSLisandro Dalcin       int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1;
1553*066ea43fSLisandro Dalcin       int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1;
1554*066ea43fSLisandro Dalcin       isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE;
1555*066ea43fSLisandro Dalcin       isHybrid  = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE;
1556*066ea43fSLisandro Dalcin       hasTetra  = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE;
1557*066ea43fSLisandro Dalcin     }
1558698a79b9SLisandro Dalcin   }
1559698a79b9SLisandro Dalcin 
1560698a79b9SLisandro Dalcin   if (parentviewer) {
1561698a79b9SLisandro Dalcin     ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr);
1562698a79b9SLisandro Dalcin   }
1563698a79b9SLisandro Dalcin 
1564*066ea43fSLisandro Dalcin   {
1565*066ea43fSLisandro Dalcin     int buf[6];
1566698a79b9SLisandro Dalcin 
1567*066ea43fSLisandro Dalcin     buf[0] = (int)dim;
1568*066ea43fSLisandro Dalcin     buf[1] = (int)order;
1569*066ea43fSLisandro Dalcin     buf[2] = periodic;
1570*066ea43fSLisandro Dalcin     buf[3] = isSimplex;
1571*066ea43fSLisandro Dalcin     buf[4] = isHybrid;
1572*066ea43fSLisandro Dalcin     buf[5] = hasTetra;
1573*066ea43fSLisandro Dalcin 
1574*066ea43fSLisandro Dalcin     ierr = MPI_Bcast(buf, 6, MPI_INT, 0, comm);CHKERRQ(ierr);
1575*066ea43fSLisandro Dalcin 
1576*066ea43fSLisandro Dalcin     dim       = buf[0];
1577*066ea43fSLisandro Dalcin     order     = buf[1];
1578*066ea43fSLisandro Dalcin     periodic  = buf[2] ? PETSC_TRUE : PETSC_FALSE;
1579*066ea43fSLisandro Dalcin     isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE;
1580*066ea43fSLisandro Dalcin     isHybrid  = buf[4] ? PETSC_TRUE : PETSC_FALSE;
1581*066ea43fSLisandro Dalcin     hasTetra  = buf[5] ? PETSC_TRUE : PETSC_FALSE;
1582dea2a3c7SStefano Zampini   }
1583dea2a3c7SStefano Zampini 
1584*066ea43fSLisandro Dalcin   if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE;
1585*066ea43fSLisandro Dalcin   if (highOrder && isHybrid) SETERRQ(comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet");
1586*066ea43fSLisandro Dalcin 
15870598e1a2SLisandro Dalcin   /* We do not want this label automatically computed, instead we fill it here */
15880598e1a2SLisandro Dalcin   ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr);
158911c56cb0SLisandro Dalcin 
1590a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
15910598e1a2SLisandro Dalcin   ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr);
15920598e1a2SLisandro Dalcin   for (cell = 0; cell < numCells; ++cell) {
15930598e1a2SLisandro Dalcin     GmshElement *elem = mesh->elements + cell;
15940598e1a2SLisandro Dalcin     DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType);
15950598e1a2SLisandro Dalcin     if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR;
15960598e1a2SLisandro Dalcin     ierr = DMPlexSetConeSize(*dm, cell, elem->numVerts);CHKERRQ(ierr);
15970598e1a2SLisandro Dalcin     ierr = DMPlexSetCellType(*dm, cell, ctype);CHKERRQ(ierr);
1598db397164SMichael Lange   }
15990598e1a2SLisandro Dalcin   for (v = numCells; v < numCells+numVerts; ++v) {
160096ca5757SLisandro Dalcin     ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr);
160196ca5757SLisandro Dalcin   }
1602331830f3SMatthew G. Knepley   ierr = DMSetUp(*dm);CHKERRQ(ierr);
160396ca5757SLisandro Dalcin 
1604a4bb7517SMichael Lange   /* Add cell-vertex connections */
16050598e1a2SLisandro Dalcin   for (cell = 0; cell < numCells; ++cell) {
16060598e1a2SLisandro Dalcin     GmshElement *elem = mesh->elements + cell;
16070598e1a2SLisandro Dalcin     for (v = 0; v < elem->numVerts; ++v) {
16080598e1a2SLisandro Dalcin       const PetscInt nn = elem->nodes[v];
16090598e1a2SLisandro Dalcin       const PetscInt vv = mesh->vertexMap[nn];
16100598e1a2SLisandro Dalcin       cone[v] = numCells + vv;
1611db397164SMichael Lange     }
16120598e1a2SLisandro Dalcin     ierr = DMPlexReorderCell(*dm, cell, cone);CHKERRQ(ierr);
16130598e1a2SLisandro Dalcin     ierr = DMPlexSetCone(*dm, cell, cone);CHKERRQ(ierr);
1614a4bb7517SMichael Lange   }
161596ca5757SLisandro Dalcin 
1616c73cfb54SMatthew G. Knepley   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
1617331830f3SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1618331830f3SMatthew G. Knepley   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1619331830f3SMatthew G. Knepley   if (interpolate) {
16205fd9971aSMatthew G. Knepley     DM idm;
1621331830f3SMatthew G. Knepley 
1622331830f3SMatthew G. Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1623331830f3SMatthew G. Knepley     ierr = DMDestroy(dm);CHKERRQ(ierr);
1624331830f3SMatthew G. Knepley     *dm  = idm;
1625331830f3SMatthew G. Knepley   }
1626917266a4SLisandro Dalcin 
1627f6c8a31fSLisandro Dalcin   if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
1628917266a4SLisandro Dalcin   if (!rank && usemarker) {
1629d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
1630d3f73514SStefano Zampini 
1631d3f73514SStefano Zampini     ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);
1632d3f73514SStefano Zampini     ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
1633d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
1634d3f73514SStefano Zampini       PetscInt suppSize;
1635d3f73514SStefano Zampini 
1636d3f73514SStefano Zampini       ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr);
1637d3f73514SStefano Zampini       if (suppSize == 1) {
1638d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
1639d3f73514SStefano Zampini 
1640d3f73514SStefano Zampini         ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
1641d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
1642d3f73514SStefano Zampini           ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr);
1643d3f73514SStefano Zampini         }
1644d3f73514SStefano Zampini         ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr);
1645d3f73514SStefano Zampini       }
1646d3f73514SStefano Zampini     }
1647d3f73514SStefano Zampini   }
164816ad7e67SMichael Lange 
164916ad7e67SMichael Lange   if (!rank) {
165011c56cb0SLisandro Dalcin     PetscInt vStart, vEnd;
1651d1a54cefSMatthew G. Knepley 
165216ad7e67SMichael Lange     ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
16530598e1a2SLisandro Dalcin     for (cell = 0, e = 0; e < numElems; ++e) {
16540598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + e;
16556e1dd89aSLawrence Mitchell 
16566e1dd89aSLawrence Mitchell       /* Create cell sets */
16570598e1a2SLisandro Dalcin       if (elem->dim == dim && dim > 0) {
16580598e1a2SLisandro Dalcin         if (elem->numTags > 0) {
16590598e1a2SLisandro Dalcin           ierr = DMSetLabelValue(*dm, "Cell Sets", cell, elem->tags[0]);CHKERRQ(ierr);
16606e1dd89aSLawrence Mitchell         }
1661917266a4SLisandro Dalcin         cell++;
16626e1dd89aSLawrence Mitchell       }
16630c070f12SSander Arens 
16640598e1a2SLisandro Dalcin       /* Create face sets */
16650598e1a2SLisandro Dalcin       if (interpolate && elem->dim == dim-1) {
16660598e1a2SLisandro Dalcin         PetscInt        joinSize;
16670598e1a2SLisandro Dalcin         const PetscInt *join = NULL;
16680598e1a2SLisandro Dalcin         /* Find the relevant facet with vertex joins */
16690598e1a2SLisandro Dalcin         for (v = 0; v < elem->numVerts; ++v) {
16700598e1a2SLisandro Dalcin           const PetscInt nn = elem->nodes[v];
16710598e1a2SLisandro Dalcin           const PetscInt vv = mesh->vertexMap[nn];
16720598e1a2SLisandro Dalcin           cone[v] = vStart + vv;
16730598e1a2SLisandro Dalcin         }
16740598e1a2SLisandro Dalcin         ierr = DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr);
16750598e1a2SLisandro Dalcin         if (joinSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Could not determine Plex facet for Gmsh element %D (Plex cell %D)", elem->id, e);
16760598e1a2SLisandro Dalcin         ierr = DMSetLabelValue(*dm, "Face Sets", join[0], elem->tags[0]);CHKERRQ(ierr);
16770598e1a2SLisandro Dalcin         ierr = DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr);
16780598e1a2SLisandro Dalcin       }
16790598e1a2SLisandro Dalcin 
16800c070f12SSander Arens       /* Create vertex sets */
16810598e1a2SLisandro Dalcin       if (elem->dim == 0) {
16820598e1a2SLisandro Dalcin         if (elem->numTags > 0) {
16830598e1a2SLisandro Dalcin           const PetscInt nn = elem->nodes[0];
16840598e1a2SLisandro Dalcin           const PetscInt vv = mesh->vertexMap[nn];
16850598e1a2SLisandro Dalcin           ierr = DMSetLabelValue(*dm, "Vertex Sets", vStart + vv, elem->tags[0]);CHKERRQ(ierr);
16860598e1a2SLisandro Dalcin         }
16870598e1a2SLisandro Dalcin       }
16880598e1a2SLisandro Dalcin     }
16890598e1a2SLisandro Dalcin   }
16900598e1a2SLisandro Dalcin 
16910598e1a2SLisandro Dalcin   if (periodic) {
16920598e1a2SLisandro Dalcin     ierr = PetscBTCreate(numVerts, &periodicVerts);CHKERRQ(ierr);
16930598e1a2SLisandro Dalcin     for (n = 0; n < numNodes; ++n) {
16940598e1a2SLisandro Dalcin       if (mesh->vertexMap[n] >= 0) {
16950598e1a2SLisandro Dalcin         if (PetscUnlikely(mesh->periodMap[n] != n)) {
16960598e1a2SLisandro Dalcin           PetscInt m = mesh->periodMap[n];
16970598e1a2SLisandro Dalcin           ierr= PetscBTSet(periodicVerts, mesh->vertexMap[n]);CHKERRQ(ierr);
16980598e1a2SLisandro Dalcin           ierr= PetscBTSet(periodicVerts, mesh->vertexMap[m]);CHKERRQ(ierr);
16990598e1a2SLisandro Dalcin         }
17000598e1a2SLisandro Dalcin       }
17010598e1a2SLisandro Dalcin     }
17020598e1a2SLisandro Dalcin     ierr = PetscBTCreate(numCells, &periodicCells);CHKERRQ(ierr);
17030598e1a2SLisandro Dalcin     for (cell = 0; cell < numCells; ++cell) {
17040598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + cell;
17050598e1a2SLisandro Dalcin       for (v = 0; v < elem->numVerts; ++v) {
17060598e1a2SLisandro Dalcin         PetscInt nn = elem->nodes[v];
17070598e1a2SLisandro Dalcin         PetscInt vv = mesh->vertexMap[nn];
17080598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) {
17090598e1a2SLisandro Dalcin           ierr = PetscBTSet(periodicCells, cell);CHKERRQ(ierr); break;
17100c070f12SSander Arens         }
17110c070f12SSander Arens       }
17120c070f12SSander Arens     }
171316ad7e67SMichael Lange   }
171416ad7e67SMichael Lange 
1715*066ea43fSLisandro Dalcin   /* Setup coordinate DM */
17160598e1a2SLisandro Dalcin   if (coordDim < 0) coordDim = dim;
17170598e1a2SLisandro Dalcin   ierr = DMSetCoordinateDim(*dm, coordDim);CHKERRQ(ierr);
1718*066ea43fSLisandro Dalcin   ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr);
1719*066ea43fSLisandro Dalcin   if (highOrder) {
1720*066ea43fSLisandro Dalcin     PetscFE         fe;
1721*066ea43fSLisandro Dalcin     PetscBool       continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
1722*066ea43fSLisandro Dalcin     PetscDTNodeType nodeType   = PETSCDTNODES_EQUISPACED;
1723*066ea43fSLisandro Dalcin 
1724*066ea43fSLisandro Dalcin     if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */
1725*066ea43fSLisandro Dalcin 
1726*066ea43fSLisandro Dalcin     ierr = GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe);CHKERRQ(ierr);
1727*066ea43fSLisandro Dalcin     ierr = PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view");CHKERRQ(ierr);
1728*066ea43fSLisandro Dalcin     ierr = DMSetField(cdm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr);
1729*066ea43fSLisandro Dalcin     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1730*066ea43fSLisandro Dalcin     ierr = DMCreateDS(cdm);CHKERRQ(ierr);
1731*066ea43fSLisandro Dalcin   }
1732*066ea43fSLisandro Dalcin 
1733*066ea43fSLisandro Dalcin   /* Create coordinates */
1734*066ea43fSLisandro Dalcin   if (highOrder) {
1735*066ea43fSLisandro Dalcin 
1736*066ea43fSLisandro Dalcin     PetscInt     maxDof = GmshNumNodes_HEX(order)*coordDim;
1737*066ea43fSLisandro Dalcin     double       *coords = mesh ? mesh->nodelist->xyz : NULL;
1738*066ea43fSLisandro Dalcin     PetscSection section;
1739*066ea43fSLisandro Dalcin     PetscScalar  *cellCoords;
1740*066ea43fSLisandro Dalcin 
1741*066ea43fSLisandro Dalcin     ierr = DMSetLocalSection(cdm, NULL);CHKERRQ(ierr);
1742*066ea43fSLisandro Dalcin     ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
1743*066ea43fSLisandro Dalcin     ierr = PetscSectionClone(coordSection, &section);CHKERRQ(ierr);
1744*066ea43fSLisandro Dalcin     ierr = DMPlexSetClosurePermutationTensor(cdm, 0, section);CHKERRQ(ierr); /* XXX Implement DMPlexSetClosurePermutationLexicographic() */
1745*066ea43fSLisandro Dalcin 
1746*066ea43fSLisandro Dalcin     ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr);
1747*066ea43fSLisandro Dalcin     ierr = PetscMalloc1(maxDof, &cellCoords);CHKERRQ(ierr);
1748*066ea43fSLisandro Dalcin     for (cell = 0; cell < numCells; ++cell) {
1749*066ea43fSLisandro Dalcin       GmshElement *elem = mesh->elements + cell;
1750*066ea43fSLisandro Dalcin       const int *lexorder = GmshCellMap[elem->cellType].lexorder();
1751*066ea43fSLisandro Dalcin       for (n = 0; n < elem->numNodes; ++n) {
1752*066ea43fSLisandro Dalcin         const PetscInt node = elem->nodes[lexorder[n]];
1753*066ea43fSLisandro Dalcin         for (d = 0; d < coordDim; ++d)
1754*066ea43fSLisandro Dalcin           cellCoords[n*coordDim+d] = (PetscReal) coords[node*3+d];
1755*066ea43fSLisandro Dalcin       }
1756*066ea43fSLisandro Dalcin       ierr = DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES);CHKERRQ(ierr);
1757*066ea43fSLisandro Dalcin     }
1758*066ea43fSLisandro Dalcin     ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
1759*066ea43fSLisandro Dalcin     ierr = PetscFree(cellCoords);CHKERRQ(ierr);
1760*066ea43fSLisandro Dalcin 
1761*066ea43fSLisandro Dalcin   } else {
1762*066ea43fSLisandro Dalcin 
1763*066ea43fSLisandro Dalcin     PetscInt    *nodeMap;
1764*066ea43fSLisandro Dalcin     double      *coords = mesh ? mesh->nodelist->xyz : NULL;
1765*066ea43fSLisandro Dalcin     PetscScalar *pointCoords;
1766*066ea43fSLisandro Dalcin 
1767*066ea43fSLisandro Dalcin     ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr);
1768331830f3SMatthew G. Knepley     ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
17690598e1a2SLisandro Dalcin     ierr = PetscSectionSetFieldComponents(coordSection, 0, coordDim);CHKERRQ(ierr);
17700598e1a2SLisandro Dalcin     if (periodic) { /* we need to localize coordinates on cells */
17710598e1a2SLisandro Dalcin       ierr = PetscSectionSetChart(coordSection, 0, numCells+numVerts);CHKERRQ(ierr);
1772f45c9589SStefano Zampini     } else {
17730598e1a2SLisandro Dalcin       ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVerts);CHKERRQ(ierr);
1774f45c9589SStefano Zampini     }
17750598e1a2SLisandro Dalcin     if (periodic) {
17760598e1a2SLisandro Dalcin       for (cell = 0; cell < numCells; ++cell) {
17770598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) {
17780598e1a2SLisandro Dalcin           GmshElement *elem = mesh->elements + cell;
17790598e1a2SLisandro Dalcin           PetscInt dof = elem->numVerts * coordDim;
17800598e1a2SLisandro Dalcin           ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr);
17810598e1a2SLisandro Dalcin           ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr);
1782f45c9589SStefano Zampini         }
1783f45c9589SStefano Zampini       }
1784f45c9589SStefano Zampini     }
17850598e1a2SLisandro Dalcin     for (v = numCells; v < numCells+numVerts; ++v) {
17860598e1a2SLisandro Dalcin       ierr = PetscSectionSetDof(coordSection, v, coordDim);CHKERRQ(ierr);
17870598e1a2SLisandro Dalcin       ierr = PetscSectionSetFieldDof(coordSection, v, 0, coordDim);CHKERRQ(ierr);
17880598e1a2SLisandro Dalcin     }
1789331830f3SMatthew G. Knepley     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1790*066ea43fSLisandro Dalcin 
1791*066ea43fSLisandro Dalcin     ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr);
1792*066ea43fSLisandro Dalcin     ierr = VecGetArray(coordinates, &pointCoords);CHKERRQ(ierr);
17930598e1a2SLisandro Dalcin     if (periodic) {
17940598e1a2SLisandro Dalcin       for (cell = 0; cell < numCells; ++cell) {
17950598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) {
17960598e1a2SLisandro Dalcin           GmshElement *elem = mesh->elements + cell;
17970598e1a2SLisandro Dalcin           PetscInt off, node;
17980598e1a2SLisandro Dalcin           for (v = 0; v < elem->numVerts; ++v)
17990598e1a2SLisandro Dalcin             cone[v] = elem->nodes[v];
1800*066ea43fSLisandro Dalcin           ierr = DMPlexReorderCell(cdm, cell, cone);CHKERRQ(ierr);
18010598e1a2SLisandro Dalcin           ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr);
18020598e1a2SLisandro Dalcin           for (v = 0; v < elem->numVerts; ++v)
18030598e1a2SLisandro Dalcin             for (node = cone[v], d = 0; d < coordDim; ++d)
1804*066ea43fSLisandro Dalcin               pointCoords[off++] = (PetscReal) coords[node*3+d];
1805331830f3SMatthew G. Knepley         }
1806331830f3SMatthew G. Knepley       }
1807331830f3SMatthew G. Knepley     }
1808*066ea43fSLisandro Dalcin     ierr = PetscMalloc1(numVerts, &nodeMap);CHKERRQ(ierr);
18090598e1a2SLisandro Dalcin     for (n = 0; n < numNodes; n++)
18100598e1a2SLisandro Dalcin       if (mesh->vertexMap[n] >= 0)
1811*066ea43fSLisandro Dalcin         nodeMap[mesh->vertexMap[n]] = n;
18120598e1a2SLisandro Dalcin     for (v = 0; v < numVerts; ++v) {
1813*066ea43fSLisandro Dalcin       PetscInt off, node = nodeMap[v];
18140598e1a2SLisandro Dalcin       ierr = PetscSectionGetOffset(coordSection, numCells + v, &off);CHKERRQ(ierr);
18150598e1a2SLisandro Dalcin       for (d = 0; d < coordDim; ++d)
1816*066ea43fSLisandro Dalcin         pointCoords[off+d] = (PetscReal) coords[node*3+d];
18170598e1a2SLisandro Dalcin     }
1818*066ea43fSLisandro Dalcin     ierr = PetscFree(nodeMap);CHKERRQ(ierr);
1819*066ea43fSLisandro Dalcin     ierr = VecRestoreArray(coordinates, &pointCoords);CHKERRQ(ierr);
1820*066ea43fSLisandro Dalcin 
1821*066ea43fSLisandro Dalcin   }
1822*066ea43fSLisandro Dalcin 
1823*066ea43fSLisandro Dalcin   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1824*066ea43fSLisandro Dalcin   ierr = VecSetBlockSize(coordinates, coordDim);CHKERRQ(ierr);
1825331830f3SMatthew G. Knepley   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
18260598e1a2SLisandro Dalcin   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
182711c56cb0SLisandro Dalcin   ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr);
182811c56cb0SLisandro Dalcin 
18290598e1a2SLisandro Dalcin   ierr = GmshMeshDestroy(&mesh);CHKERRQ(ierr);
18300598e1a2SLisandro Dalcin   ierr = PetscBTDestroy(&periodicVerts);CHKERRQ(ierr);
18310598e1a2SLisandro Dalcin   ierr = PetscBTDestroy(&periodicCells);CHKERRQ(ierr);
183211c56cb0SLisandro Dalcin 
1833*066ea43fSLisandro Dalcin   if (highOrder && project)  {
1834*066ea43fSLisandro Dalcin     PetscFE         fe;
1835*066ea43fSLisandro Dalcin     const char      prefix[]   = "dm_plex_gmsh_project_";
1836*066ea43fSLisandro Dalcin     PetscBool       continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
1837*066ea43fSLisandro Dalcin     PetscDTNodeType nodeType   = PETSCDTNODES_GAUSSJACOBI;
1838*066ea43fSLisandro Dalcin 
1839*066ea43fSLisandro Dalcin     if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */
1840*066ea43fSLisandro Dalcin 
1841*066ea43fSLisandro Dalcin     ierr = GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe);CHKERRQ(ierr);
1842*066ea43fSLisandro Dalcin     ierr = PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view");CHKERRQ(ierr);
1843*066ea43fSLisandro Dalcin     ierr = DMProjectCoordinates(*dm, fe);CHKERRQ(ierr);
1844*066ea43fSLisandro Dalcin     ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
1845*066ea43fSLisandro Dalcin   }
1846*066ea43fSLisandro Dalcin 
18470598e1a2SLisandro Dalcin   ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr);
1848331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
1849331830f3SMatthew G. Knepley }
1850