xref: /petsc/src/dm/impls/plex/plexgmsh.c (revision 5f5cd6d5074be84d275388b323090bd5318d3f81)
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 
5066ea43fSLisandro Dalcin #include <../src/dm/impls/plex/gmshlex.h>
6066ea43fSLisandro Dalcin 
7066ea43fSLisandro Dalcin #define GMSH_LEXORDER_ITEM(T, p)                                   \
8066ea43fSLisandro Dalcin static int *GmshLexOrder_##T##_##p(void)                           \
9066ea43fSLisandro Dalcin {                                                                  \
10066ea43fSLisandro Dalcin   static int Gmsh_LexOrder_##T##_##p[GmshNumNodes_##T(p)] = {-1};  \
11066ea43fSLisandro Dalcin   int *lex = Gmsh_LexOrder_##T##_##p;                              \
12066ea43fSLisandro Dalcin   if (lex[0] == -1) (void)GmshLexOrder_##T(p, lex, 0);             \
13066ea43fSLisandro Dalcin   return lex;                                                      \
14066ea43fSLisandro Dalcin }
15066ea43fSLisandro Dalcin 
16b9bf55e5SMatthew G. Knepley static int *GmshLexOrder_QUA_2_Serendipity(void)
17b9bf55e5SMatthew G. Knepley {
18b9bf55e5SMatthew G. Knepley   static int Gmsh_LexOrder_QUA_2_Serendipity[9] = {-1};
19b9bf55e5SMatthew G. Knepley   int *lex = Gmsh_LexOrder_QUA_2_Serendipity;
20b9bf55e5SMatthew G. Knepley   if (lex[0] == -1) {
21b9bf55e5SMatthew G. Knepley     /* Vertices */
22b9bf55e5SMatthew G. Knepley     lex[0] = 0; lex[2] = 1; lex[8] = 2; lex[6] = 3;
23b9bf55e5SMatthew G. Knepley     /* Edges */
24b9bf55e5SMatthew G. Knepley     lex[1] = 4; lex[5] = 5; lex[7] = 6; lex[3] = 7;
25b9bf55e5SMatthew G. Knepley     /* Cell */
26b9bf55e5SMatthew G. Knepley     lex[4] = -8;
27b9bf55e5SMatthew G. Knepley   }
28b9bf55e5SMatthew G. Knepley   return lex;
29b9bf55e5SMatthew G. Knepley }
30b9bf55e5SMatthew G. Knepley 
31b9bf55e5SMatthew G. Knepley static int *GmshLexOrder_HEX_2_Serendipity(void)
32b9bf55e5SMatthew G. Knepley {
33b9bf55e5SMatthew G. Knepley   static int Gmsh_LexOrder_HEX_2_Serendipity[27] = {-1};
34b9bf55e5SMatthew G. Knepley   int *lex = Gmsh_LexOrder_HEX_2_Serendipity;
35b9bf55e5SMatthew G. Knepley   if (lex[0] == -1) {
36b9bf55e5SMatthew G. Knepley     /* Vertices */
37b9bf55e5SMatthew G. Knepley     lex[ 0] =   0; lex[ 2] =   1; lex[ 8] =   2; lex[ 6] =   3;
38b9bf55e5SMatthew G. Knepley     lex[18] =   4; lex[20] =   5; lex[26] =   6; lex[24] =   7;
39b9bf55e5SMatthew G. Knepley     /* Edges */
40b9bf55e5SMatthew G. Knepley     lex[ 1] =   8; lex[ 3] =   9; lex[ 9] =  10; lex[ 5] =  11;
41b9bf55e5SMatthew G. Knepley     lex[11] =  12; lex[ 7] =  13; lex[17] =  14; lex[15] =  15;
42b9bf55e5SMatthew G. Knepley     lex[19] =  16; lex[21] =  17; lex[23] =  18; lex[25] =  19;
43b9bf55e5SMatthew G. Knepley     /* Faces */
44b9bf55e5SMatthew G. Knepley     lex[ 4] = -20; lex[10] = -21; lex[12] = -22; lex[14] = -23;
45b9bf55e5SMatthew G. Knepley     lex[16] = -24; lex[22] = -25;
46b9bf55e5SMatthew G. Knepley     /* Cell */
47b9bf55e5SMatthew G. Knepley     lex[13] = -26;
48b9bf55e5SMatthew G. Knepley   }
49b9bf55e5SMatthew G. Knepley   return lex;
50b9bf55e5SMatthew G. Knepley }
51b9bf55e5SMatthew G. Knepley 
52066ea43fSLisandro Dalcin #define GMSH_LEXORDER_LIST(T) \
53066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  1)     \
54066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  2)     \
55066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  3)     \
56066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  4)     \
57066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  5)     \
58066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  6)     \
59066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  7)     \
60066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  8)     \
61066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T,  9)     \
62066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 10)
63066ea43fSLisandro Dalcin 
64066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(VTX, 0)
65066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(SEG)
66066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TRI)
67066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(QUA)
68066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TET)
69066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(HEX)
70066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PRI)
71066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PYR)
72066ea43fSLisandro Dalcin 
73066ea43fSLisandro Dalcin typedef enum {
74066ea43fSLisandro Dalcin   GMSH_VTX = 0,
75066ea43fSLisandro Dalcin   GMSH_SEG = 1,
76066ea43fSLisandro Dalcin   GMSH_TRI = 2,
77066ea43fSLisandro Dalcin   GMSH_QUA = 3,
78066ea43fSLisandro Dalcin   GMSH_TET = 4,
79066ea43fSLisandro Dalcin   GMSH_HEX = 5,
80066ea43fSLisandro Dalcin   GMSH_PRI = 6,
81066ea43fSLisandro Dalcin   GMSH_PYR = 7,
82066ea43fSLisandro Dalcin   GMSH_NUM_POLYTOPES = 8
83066ea43fSLisandro Dalcin } GmshPolytopeType;
84066ea43fSLisandro Dalcin 
85de78e4feSLisandro Dalcin typedef struct {
860598e1a2SLisandro Dalcin   int   cellType;
87066ea43fSLisandro Dalcin   int   polytope;
880598e1a2SLisandro Dalcin   int   dim;
890598e1a2SLisandro Dalcin   int   order;
90066ea43fSLisandro Dalcin   int   numVerts;
910598e1a2SLisandro Dalcin   int   numNodes;
92066ea43fSLisandro Dalcin   int* (*lexorder)(void);
930598e1a2SLisandro Dalcin } GmshCellInfo;
940598e1a2SLisandro Dalcin 
95066ea43fSLisandro Dalcin #define GmshCellEntry(cellType, polytope, dim, order) \
96066ea43fSLisandro Dalcin   {cellType, GMSH_##polytope, dim, order, \
97066ea43fSLisandro Dalcin    GmshNumNodes_##polytope(1), \
98066ea43fSLisandro Dalcin    GmshNumNodes_##polytope(order), \
99066ea43fSLisandro Dalcin    GmshLexOrder_##polytope##_##order}
1000598e1a2SLisandro Dalcin 
1010598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = {
102066ea43fSLisandro Dalcin   GmshCellEntry( 15, VTX, 0,  0),
1030598e1a2SLisandro Dalcin 
104066ea43fSLisandro Dalcin   GmshCellEntry(  1, SEG, 1,  1),
105066ea43fSLisandro Dalcin   GmshCellEntry(  8, SEG, 1,  2),
106066ea43fSLisandro Dalcin   GmshCellEntry( 26, SEG, 1,  3),
107066ea43fSLisandro Dalcin   GmshCellEntry( 27, SEG, 1,  4),
108066ea43fSLisandro Dalcin   GmshCellEntry( 28, SEG, 1,  5),
109066ea43fSLisandro Dalcin   GmshCellEntry( 62, SEG, 1,  6),
110066ea43fSLisandro Dalcin   GmshCellEntry( 63, SEG, 1,  7),
111066ea43fSLisandro Dalcin   GmshCellEntry( 64, SEG, 1,  8),
112066ea43fSLisandro Dalcin   GmshCellEntry( 65, SEG, 1,  9),
113066ea43fSLisandro Dalcin   GmshCellEntry( 66, SEG, 1, 10),
1140598e1a2SLisandro Dalcin 
115066ea43fSLisandro Dalcin   GmshCellEntry(  2, TRI, 2,  1),
116066ea43fSLisandro Dalcin   GmshCellEntry(  9, TRI, 2,  2),
117066ea43fSLisandro Dalcin   GmshCellEntry( 21, TRI, 2,  3),
118066ea43fSLisandro Dalcin   GmshCellEntry( 23, TRI, 2,  4),
119066ea43fSLisandro Dalcin   GmshCellEntry( 25, TRI, 2,  5),
120066ea43fSLisandro Dalcin   GmshCellEntry( 42, TRI, 2,  6),
121066ea43fSLisandro Dalcin   GmshCellEntry( 43, TRI, 2,  7),
122066ea43fSLisandro Dalcin   GmshCellEntry( 44, TRI, 2,  8),
123066ea43fSLisandro Dalcin   GmshCellEntry( 45, TRI, 2,  9),
124066ea43fSLisandro Dalcin   GmshCellEntry( 46, TRI, 2, 10),
1250598e1a2SLisandro Dalcin 
126066ea43fSLisandro Dalcin   GmshCellEntry(  3, QUA, 2,  1),
127066ea43fSLisandro Dalcin   GmshCellEntry( 10, QUA, 2,  2),
128b9bf55e5SMatthew G. Knepley   {16, GMSH_QUA, 2, 2, 4, 8, GmshLexOrder_QUA_2_Serendipity},
129066ea43fSLisandro Dalcin   GmshCellEntry( 36, QUA, 2,  3),
130066ea43fSLisandro Dalcin   GmshCellEntry( 37, QUA, 2,  4),
131066ea43fSLisandro Dalcin   GmshCellEntry( 38, QUA, 2,  5),
132066ea43fSLisandro Dalcin   GmshCellEntry( 47, QUA, 2,  6),
133066ea43fSLisandro Dalcin   GmshCellEntry( 48, QUA, 2,  7),
134066ea43fSLisandro Dalcin   GmshCellEntry( 49, QUA, 2,  8),
135066ea43fSLisandro Dalcin   GmshCellEntry( 50, QUA, 2,  9),
136066ea43fSLisandro Dalcin   GmshCellEntry( 51, QUA, 2, 10),
1370598e1a2SLisandro Dalcin 
138066ea43fSLisandro Dalcin   GmshCellEntry(  4, TET, 3,  1),
139066ea43fSLisandro Dalcin   GmshCellEntry( 11, TET, 3,  2),
140066ea43fSLisandro Dalcin   GmshCellEntry( 29, TET, 3,  3),
141066ea43fSLisandro Dalcin   GmshCellEntry( 30, TET, 3,  4),
142066ea43fSLisandro Dalcin   GmshCellEntry( 31, TET, 3,  5),
143066ea43fSLisandro Dalcin   GmshCellEntry( 71, TET, 3,  6),
144066ea43fSLisandro Dalcin   GmshCellEntry( 72, TET, 3,  7),
145066ea43fSLisandro Dalcin   GmshCellEntry( 73, TET, 3,  8),
146066ea43fSLisandro Dalcin   GmshCellEntry( 74, TET, 3,  9),
147066ea43fSLisandro Dalcin   GmshCellEntry( 75, TET, 3, 10),
1480598e1a2SLisandro Dalcin 
149066ea43fSLisandro Dalcin   GmshCellEntry(  5, HEX, 3,  1),
150066ea43fSLisandro Dalcin   GmshCellEntry( 12, HEX, 3,  2),
151b9bf55e5SMatthew G. Knepley   {17, GMSH_HEX, 3, 2, 8, 20, GmshLexOrder_HEX_2_Serendipity},
152066ea43fSLisandro Dalcin   GmshCellEntry( 92, HEX, 3,  3),
153066ea43fSLisandro Dalcin   GmshCellEntry( 93, HEX, 3,  4),
154066ea43fSLisandro Dalcin   GmshCellEntry( 94, HEX, 3,  5),
155066ea43fSLisandro Dalcin   GmshCellEntry( 95, HEX, 3,  6),
156066ea43fSLisandro Dalcin   GmshCellEntry( 96, HEX, 3,  7),
157066ea43fSLisandro Dalcin   GmshCellEntry( 97, HEX, 3,  8),
158066ea43fSLisandro Dalcin   GmshCellEntry( 98, HEX, 3,  9),
159066ea43fSLisandro Dalcin   GmshCellEntry( -1, HEX, 3, 10),
1600598e1a2SLisandro Dalcin 
161066ea43fSLisandro Dalcin   GmshCellEntry(  6, PRI, 3,  1),
162066ea43fSLisandro Dalcin   GmshCellEntry( 13, PRI, 3,  2),
163066ea43fSLisandro Dalcin   GmshCellEntry( 90, PRI, 3,  3),
164066ea43fSLisandro Dalcin   GmshCellEntry( 91, PRI, 3,  4),
165066ea43fSLisandro Dalcin   GmshCellEntry(106, PRI, 3,  5),
166066ea43fSLisandro Dalcin   GmshCellEntry(107, PRI, 3,  6),
167066ea43fSLisandro Dalcin   GmshCellEntry(108, PRI, 3,  7),
168066ea43fSLisandro Dalcin   GmshCellEntry(109, PRI, 3,  8),
169066ea43fSLisandro Dalcin   GmshCellEntry(110, PRI, 3,  9),
170066ea43fSLisandro Dalcin   GmshCellEntry( -1, PRI, 3, 10),
1710598e1a2SLisandro Dalcin 
172066ea43fSLisandro Dalcin   GmshCellEntry(  7, PYR, 3,  1),
173066ea43fSLisandro Dalcin   GmshCellEntry( 14, PYR, 3,  2),
174066ea43fSLisandro Dalcin   GmshCellEntry(118, PYR, 3,  3),
175066ea43fSLisandro Dalcin   GmshCellEntry(119, PYR, 3,  4),
176066ea43fSLisandro Dalcin   GmshCellEntry(120, PYR, 3,  5),
177066ea43fSLisandro Dalcin   GmshCellEntry(121, PYR, 3,  6),
178066ea43fSLisandro Dalcin   GmshCellEntry(122, PYR, 3,  7),
179066ea43fSLisandro Dalcin   GmshCellEntry(123, PYR, 3,  8),
180066ea43fSLisandro Dalcin   GmshCellEntry(124, PYR, 3,  9),
181066ea43fSLisandro Dalcin   GmshCellEntry( -1, PYR, 3, 10)
1820598e1a2SLisandro Dalcin 
1830598e1a2SLisandro Dalcin #if 0
184066ea43fSLisandro Dalcin   {20, GMSH_TRI, 2, 3, 3,  9, NULL},
185066ea43fSLisandro Dalcin   {18, GMSH_PRI, 3, 2, 6, 15, NULL},
186066ea43fSLisandro Dalcin   {19, GMSH_PYR, 3, 2, 5, 13, NULL},
1870598e1a2SLisandro Dalcin #endif
1880598e1a2SLisandro Dalcin };
1890598e1a2SLisandro Dalcin 
1900598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150];
1910598e1a2SLisandro Dalcin 
1920598e1a2SLisandro Dalcin static PetscErrorCode GmshCellInfoSetUp(void)
1930598e1a2SLisandro Dalcin {
1940598e1a2SLisandro Dalcin   size_t           i, n;
1950598e1a2SLisandro Dalcin   static PetscBool called = PETSC_FALSE;
1960598e1a2SLisandro Dalcin 
1970598e1a2SLisandro Dalcin   if (called) return 0;
1980598e1a2SLisandro Dalcin   PetscFunctionBegin;
1990598e1a2SLisandro Dalcin   called = PETSC_TRUE;
200dd39110bSPierre Jolivet   n = PETSC_STATIC_ARRAY_LENGTH(GmshCellMap);
2010598e1a2SLisandro Dalcin   for (i = 0; i < n; ++i) {
2020598e1a2SLisandro Dalcin     GmshCellMap[i].cellType = -1;
203066ea43fSLisandro Dalcin     GmshCellMap[i].polytope = -1;
2040598e1a2SLisandro Dalcin   }
205dd39110bSPierre Jolivet   n = PETSC_STATIC_ARRAY_LENGTH(GmshCellTable);
206066ea43fSLisandro Dalcin   for (i = 0; i < n; ++i) {
207066ea43fSLisandro Dalcin     if (GmshCellTable[i].cellType <= 0) continue;
208066ea43fSLisandro Dalcin     GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i];
209066ea43fSLisandro Dalcin   }
2100598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
2110598e1a2SLisandro Dalcin }
2120598e1a2SLisandro Dalcin 
21394bad497SJacob Faibussowitsch #define GmshCellTypeCheck(ct) PetscMacroReturnStandard(                                        \
2140598e1a2SLisandro Dalcin     const int _ct_ = (int)ct;                                                                  \
215dd39110bSPierre Jolivet     PetscCheck(_ct_ >= 0 && _ct_ < (int)PETSC_STATIC_ARRAY_LENGTH(GmshCellMap), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); \
21694bad497SJacob Faibussowitsch     PetscCheck(GmshCellMap[_ct_].cellType == _ct_, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \
21794bad497SJacob Faibussowitsch     PetscCheck(GmshCellMap[_ct_].polytope != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \
21894bad497SJacob Faibussowitsch   )
2190598e1a2SLisandro Dalcin 
2200598e1a2SLisandro Dalcin typedef struct {
221698a79b9SLisandro Dalcin   PetscViewer  viewer;
222698a79b9SLisandro Dalcin   int          fileFormat;
223698a79b9SLisandro Dalcin   int          dataSize;
224698a79b9SLisandro Dalcin   PetscBool    binary;
225698a79b9SLisandro Dalcin   PetscBool    byteSwap;
226698a79b9SLisandro Dalcin   size_t       wlen;
227698a79b9SLisandro Dalcin   void        *wbuf;
228698a79b9SLisandro Dalcin   size_t       slen;
229698a79b9SLisandro Dalcin   void        *sbuf;
2300598e1a2SLisandro Dalcin   PetscInt    *nbuf;
2310598e1a2SLisandro Dalcin   PetscInt     nodeStart;
2320598e1a2SLisandro Dalcin   PetscInt     nodeEnd;
23333a088b6SMatthew G. Knepley   PetscInt    *nodeMap;
234698a79b9SLisandro Dalcin } GmshFile;
235de78e4feSLisandro Dalcin 
236698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf)
237de78e4feSLisandro Dalcin {
238de78e4feSLisandro Dalcin   size_t         size = count * eltsize;
23911c56cb0SLisandro Dalcin 
24011c56cb0SLisandro Dalcin   PetscFunctionBegin;
241698a79b9SLisandro Dalcin   if (gmsh->wlen < size) {
2429566063dSJacob Faibussowitsch     PetscCall(PetscFree(gmsh->wbuf));
2439566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(size, &gmsh->wbuf));
244698a79b9SLisandro Dalcin     gmsh->wlen = size;
245de78e4feSLisandro Dalcin   }
246698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->wbuf : NULL;
247de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
248de78e4feSLisandro Dalcin }
249de78e4feSLisandro Dalcin 
250698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf)
251698a79b9SLisandro Dalcin {
252698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
253698a79b9SLisandro Dalcin   size_t         size = count * dataSize;
254698a79b9SLisandro Dalcin 
255698a79b9SLisandro Dalcin   PetscFunctionBegin;
256698a79b9SLisandro Dalcin   if (gmsh->slen < size) {
2579566063dSJacob Faibussowitsch     PetscCall(PetscFree(gmsh->sbuf));
2589566063dSJacob Faibussowitsch     PetscCall(PetscMalloc(size, &gmsh->sbuf));
259698a79b9SLisandro Dalcin     gmsh->slen = size;
260698a79b9SLisandro Dalcin   }
261698a79b9SLisandro Dalcin   *(void**)buf = size ? gmsh->sbuf : NULL;
262698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
263698a79b9SLisandro Dalcin }
264698a79b9SLisandro Dalcin 
265698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype)
266de78e4feSLisandro Dalcin {
267de78e4feSLisandro Dalcin   PetscFunctionBegin;
2689566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype));
2699566063dSJacob Faibussowitsch   if (gmsh->byteSwap) PetscCall(PetscByteSwap(buf, dtype, count));
270698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
271698a79b9SLisandro Dalcin }
272698a79b9SLisandro Dalcin 
273698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count)
274698a79b9SLisandro Dalcin {
275698a79b9SLisandro Dalcin   PetscFunctionBegin;
2769566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING));
277698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
278698a79b9SLisandro Dalcin }
279698a79b9SLisandro Dalcin 
280d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match)
281d6689ca9SLisandro Dalcin {
282d6689ca9SLisandro Dalcin   PetscFunctionBegin;
2839566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(line, Section, match));
284d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
285d6689ca9SLisandro Dalcin }
286d6689ca9SLisandro Dalcin 
287d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN])
288d6689ca9SLisandro Dalcin {
289d6689ca9SLisandro Dalcin   PetscBool      match;
290d6689ca9SLisandro Dalcin 
291d6689ca9SLisandro Dalcin   PetscFunctionBegin;
2929566063dSJacob Faibussowitsch   PetscCall(GmshMatch(gmsh, Section, line, &match));
29328b400f6SJacob Faibussowitsch   PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s",Section);
294d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
295d6689ca9SLisandro Dalcin }
296d6689ca9SLisandro Dalcin 
297d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN])
298d6689ca9SLisandro Dalcin {
299d6689ca9SLisandro Dalcin   PetscBool      match;
300d6689ca9SLisandro Dalcin 
301d6689ca9SLisandro Dalcin   PetscFunctionBegin;
302d6689ca9SLisandro Dalcin   while (PETSC_TRUE) {
3039566063dSJacob Faibussowitsch     PetscCall(GmshReadString(gmsh, line, 1));
3049566063dSJacob Faibussowitsch     PetscCall(GmshMatch(gmsh, "$Comments", line, &match));
305d6689ca9SLisandro Dalcin     if (!match) break;
306d6689ca9SLisandro Dalcin     while (PETSC_TRUE) {
3079566063dSJacob Faibussowitsch       PetscCall(GmshReadString(gmsh, line, 1));
3089566063dSJacob Faibussowitsch       PetscCall(GmshMatch(gmsh, "$EndComments", line, &match));
309d6689ca9SLisandro Dalcin       if (match) break;
310d6689ca9SLisandro Dalcin     }
311d6689ca9SLisandro Dalcin   }
312d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
313d6689ca9SLisandro Dalcin }
314d6689ca9SLisandro Dalcin 
315d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN])
316d6689ca9SLisandro Dalcin {
317d6689ca9SLisandro Dalcin   PetscFunctionBegin;
3189566063dSJacob Faibussowitsch   PetscCall(GmshReadString(gmsh, line, 1));
3199566063dSJacob Faibussowitsch   PetscCall(GmshExpect(gmsh, EndSection, line));
320d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
321d6689ca9SLisandro Dalcin }
322d6689ca9SLisandro Dalcin 
323698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count)
324698a79b9SLisandro Dalcin {
325698a79b9SLisandro Dalcin   PetscInt       i;
326698a79b9SLisandro Dalcin   size_t         dataSize = (size_t)gmsh->dataSize;
327698a79b9SLisandro Dalcin 
328698a79b9SLisandro Dalcin   PetscFunctionBegin;
329698a79b9SLisandro Dalcin   if (dataSize == sizeof(PetscInt)) {
3309566063dSJacob Faibussowitsch     PetscCall(GmshRead(gmsh, buf, count, PETSC_INT));
331698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(int)) {
332698a79b9SLisandro Dalcin     int *ibuf = NULL;
3339566063dSJacob Faibussowitsch     PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
3349566063dSJacob Faibussowitsch     PetscCall(GmshRead(gmsh, ibuf, count, PETSC_ENUM));
335698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
336698a79b9SLisandro Dalcin   } else  if (dataSize == sizeof(long)) {
337698a79b9SLisandro Dalcin     long *ibuf = NULL;
3389566063dSJacob Faibussowitsch     PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
3399566063dSJacob Faibussowitsch     PetscCall(GmshRead(gmsh, ibuf, count, PETSC_LONG));
340698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
341698a79b9SLisandro Dalcin   } else if (dataSize == sizeof(PetscInt64)) {
342698a79b9SLisandro Dalcin     PetscInt64 *ibuf = NULL;
3439566063dSJacob Faibussowitsch     PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
3449566063dSJacob Faibussowitsch     PetscCall(GmshRead(gmsh, ibuf, count, PETSC_INT64));
345698a79b9SLisandro Dalcin     for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
346698a79b9SLisandro Dalcin   }
347698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
348698a79b9SLisandro Dalcin }
349698a79b9SLisandro Dalcin 
350698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count)
351698a79b9SLisandro Dalcin {
352698a79b9SLisandro Dalcin   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(GmshRead(gmsh, buf, count, PETSC_ENUM));
354698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
355698a79b9SLisandro Dalcin }
356698a79b9SLisandro Dalcin 
357698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count)
358698a79b9SLisandro Dalcin {
359698a79b9SLisandro Dalcin   PetscFunctionBegin;
3609566063dSJacob Faibussowitsch   PetscCall(GmshRead(gmsh, buf, count, PETSC_DOUBLE));
361de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
362de78e4feSLisandro Dalcin }
363de78e4feSLisandro Dalcin 
364de78e4feSLisandro Dalcin typedef struct {
3650598e1a2SLisandro Dalcin   PetscInt id;       /* Entity ID */
3660598e1a2SLisandro Dalcin   PetscInt dim;      /* Dimension */
367de78e4feSLisandro Dalcin   double   bbox[6];  /* Bounding box */
368de78e4feSLisandro Dalcin   PetscInt numTags;  /* Size of tag array */
369de78e4feSLisandro Dalcin   int      tags[4];  /* Tag array */
370de78e4feSLisandro Dalcin } GmshEntity;
371de78e4feSLisandro Dalcin 
372de78e4feSLisandro Dalcin typedef struct {
373730356f6SLisandro Dalcin   GmshEntity *entity[4];
374730356f6SLisandro Dalcin   PetscHMapI  entityMap[4];
375730356f6SLisandro Dalcin } GmshEntities;
376730356f6SLisandro Dalcin 
377698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities)
378730356f6SLisandro Dalcin {
379730356f6SLisandro Dalcin   PetscInt       dim;
380730356f6SLisandro Dalcin 
381730356f6SLisandro Dalcin   PetscFunctionBegin;
3829566063dSJacob Faibussowitsch   PetscCall(PetscNew(entities));
383730356f6SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
3849566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(count[dim], &(*entities)->entity[dim]));
3859566063dSJacob Faibussowitsch     PetscCall(PetscHMapICreate(&(*entities)->entityMap[dim]));
386730356f6SLisandro Dalcin   }
387730356f6SLisandro Dalcin   PetscFunctionReturn(0);
388730356f6SLisandro Dalcin }
389730356f6SLisandro Dalcin 
3900598e1a2SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities)
3910598e1a2SLisandro Dalcin {
3920598e1a2SLisandro Dalcin   PetscInt       dim;
3930598e1a2SLisandro Dalcin 
3940598e1a2SLisandro Dalcin   PetscFunctionBegin;
3950598e1a2SLisandro Dalcin   if (!*entities) PetscFunctionReturn(0);
3960598e1a2SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
3979566063dSJacob Faibussowitsch     PetscCall(PetscFree((*entities)->entity[dim]));
3989566063dSJacob Faibussowitsch     PetscCall(PetscHMapIDestroy(&(*entities)->entityMap[dim]));
3990598e1a2SLisandro Dalcin   }
4009566063dSJacob Faibussowitsch   PetscCall(PetscFree((*entities)));
4010598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4020598e1a2SLisandro Dalcin }
4030598e1a2SLisandro Dalcin 
404730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity)
405730356f6SLisandro Dalcin {
406730356f6SLisandro Dalcin   PetscFunctionBegin;
4079566063dSJacob Faibussowitsch   PetscCall(PetscHMapISet(entities->entityMap[dim], eid, index));
408730356f6SLisandro Dalcin   entities->entity[dim][index].dim = dim;
409730356f6SLisandro Dalcin   entities->entity[dim][index].id  = eid;
410730356f6SLisandro Dalcin   if (entity) *entity = &entities->entity[dim][index];
411730356f6SLisandro Dalcin   PetscFunctionReturn(0);
412730356f6SLisandro Dalcin }
413730356f6SLisandro Dalcin 
414730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity)
415730356f6SLisandro Dalcin {
416730356f6SLisandro Dalcin   PetscInt       index;
417730356f6SLisandro Dalcin 
418730356f6SLisandro Dalcin   PetscFunctionBegin;
4199566063dSJacob Faibussowitsch   PetscCall(PetscHMapIGet(entities->entityMap[dim], eid, &index));
420730356f6SLisandro Dalcin   *entity = &entities->entity[dim][index];
421730356f6SLisandro Dalcin   PetscFunctionReturn(0);
422730356f6SLisandro Dalcin }
423730356f6SLisandro Dalcin 
4240598e1a2SLisandro Dalcin typedef struct {
4250598e1a2SLisandro Dalcin   PetscInt *id;  /* Node IDs */
4260598e1a2SLisandro Dalcin   double   *xyz; /* Coordinates */
42781a1af93SMatthew G. Knepley   PetscInt *tag; /* Physical tag */
4280598e1a2SLisandro Dalcin } GmshNodes;
4290598e1a2SLisandro Dalcin 
4300598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes)
431730356f6SLisandro Dalcin {
432730356f6SLisandro Dalcin   PetscFunctionBegin;
4339566063dSJacob Faibussowitsch   PetscCall(PetscNew(nodes));
4349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(count*1, &(*nodes)->id));
4359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(count*3, &(*nodes)->xyz));
4369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(count*1, &(*nodes)->tag));
4370598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
438730356f6SLisandro Dalcin }
4390598e1a2SLisandro Dalcin 
4400598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes)
4410598e1a2SLisandro Dalcin {
4420598e1a2SLisandro Dalcin   PetscFunctionBegin;
4430598e1a2SLisandro Dalcin   if (!*nodes) PetscFunctionReturn(0);
4449566063dSJacob Faibussowitsch   PetscCall(PetscFree((*nodes)->id));
4459566063dSJacob Faibussowitsch   PetscCall(PetscFree((*nodes)->xyz));
4469566063dSJacob Faibussowitsch   PetscCall(PetscFree((*nodes)->tag));
4479566063dSJacob Faibussowitsch   PetscCall(PetscFree((*nodes)));
448730356f6SLisandro Dalcin   PetscFunctionReturn(0);
449730356f6SLisandro Dalcin }
450730356f6SLisandro Dalcin 
451730356f6SLisandro Dalcin typedef struct {
4520598e1a2SLisandro Dalcin   PetscInt id;       /* Element ID */
4530598e1a2SLisandro Dalcin   PetscInt dim;      /* Dimension */
454de78e4feSLisandro Dalcin   PetscInt cellType; /* Cell type */
4550598e1a2SLisandro Dalcin   PetscInt numVerts; /* Size of vertex array */
456de78e4feSLisandro Dalcin   PetscInt numNodes; /* Size of node array */
4570598e1a2SLisandro Dalcin   PetscInt *nodes;   /* Vertex/Node array */
45881a1af93SMatthew G. Knepley   PetscInt numTags;  /* Size of physical tag array */
45981a1af93SMatthew G. Knepley   int      tags[4];  /* Physical tag array */
460de78e4feSLisandro Dalcin } GmshElement;
461de78e4feSLisandro Dalcin 
4620598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements)
4630598e1a2SLisandro Dalcin {
4640598e1a2SLisandro Dalcin   PetscFunctionBegin;
4659566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(count, elements));
4660598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4670598e1a2SLisandro Dalcin }
4680598e1a2SLisandro Dalcin 
4690598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsDestroy(GmshElement **elements)
4700598e1a2SLisandro Dalcin {
4710598e1a2SLisandro Dalcin   PetscFunctionBegin;
4720598e1a2SLisandro Dalcin   if (!*elements) PetscFunctionReturn(0);
4739566063dSJacob Faibussowitsch   PetscCall(PetscFree(*elements));
4740598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
4750598e1a2SLisandro Dalcin }
4760598e1a2SLisandro Dalcin 
4770598e1a2SLisandro Dalcin typedef struct {
4780598e1a2SLisandro Dalcin   PetscInt       dim;
479066ea43fSLisandro Dalcin   PetscInt       order;
4800598e1a2SLisandro Dalcin   GmshEntities  *entities;
4810598e1a2SLisandro Dalcin   PetscInt       numNodes;
4820598e1a2SLisandro Dalcin   GmshNodes     *nodelist;
4830598e1a2SLisandro Dalcin   PetscInt       numElems;
4840598e1a2SLisandro Dalcin   GmshElement   *elements;
4850598e1a2SLisandro Dalcin   PetscInt       numVerts;
4860598e1a2SLisandro Dalcin   PetscInt       numCells;
4870598e1a2SLisandro Dalcin   PetscInt      *periodMap;
4880598e1a2SLisandro Dalcin   PetscInt      *vertexMap;
4890598e1a2SLisandro Dalcin   PetscSegBuffer segbuf;
490a45dabc8SMatthew G. Knepley   PetscInt       numRegions;
491a45dabc8SMatthew G. Knepley   PetscInt      *regionTags;
492a45dabc8SMatthew G. Knepley   char         **regionNames;
4930598e1a2SLisandro Dalcin } GmshMesh;
4940598e1a2SLisandro Dalcin 
4950598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshCreate(GmshMesh **mesh)
4960598e1a2SLisandro Dalcin {
4970598e1a2SLisandro Dalcin   PetscFunctionBegin;
4989566063dSJacob Faibussowitsch   PetscCall(PetscNew(mesh));
4999566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf));
5000598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
5010598e1a2SLisandro Dalcin }
5020598e1a2SLisandro Dalcin 
5030598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh)
5040598e1a2SLisandro Dalcin {
505a45dabc8SMatthew G. Knepley   PetscInt       r;
5060598e1a2SLisandro Dalcin 
5070598e1a2SLisandro Dalcin   PetscFunctionBegin;
5080598e1a2SLisandro Dalcin   if (!*mesh) PetscFunctionReturn(0);
5099566063dSJacob Faibussowitsch   PetscCall(GmshEntitiesDestroy(&(*mesh)->entities));
5109566063dSJacob Faibussowitsch   PetscCall(GmshNodesDestroy(&(*mesh)->nodelist));
5119566063dSJacob Faibussowitsch   PetscCall(GmshElementsDestroy(&(*mesh)->elements));
5129566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mesh)->periodMap));
5139566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mesh)->vertexMap));
5149566063dSJacob Faibussowitsch   PetscCall(PetscSegBufferDestroy(&(*mesh)->segbuf));
5159566063dSJacob Faibussowitsch   for (r = 0; r < (*mesh)->numRegions; ++r) PetscCall(PetscFree((*mesh)->regionNames[r]));
5169566063dSJacob Faibussowitsch   PetscCall(PetscFree2((*mesh)->regionTags, (*mesh)->regionNames));
5179566063dSJacob Faibussowitsch   PetscCall(PetscFree((*mesh)));
5180598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
5190598e1a2SLisandro Dalcin }
5200598e1a2SLisandro Dalcin 
5210598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh)
522de78e4feSLisandro Dalcin {
523698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
524698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
525de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
5260598e1a2SLisandro Dalcin   int            n, num, nid, snum;
5270598e1a2SLisandro Dalcin   GmshNodes      *nodes;
528de78e4feSLisandro Dalcin 
529de78e4feSLisandro Dalcin   PetscFunctionBegin;
5309566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
531698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
53208401ef6SPierre Jolivet   PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
5339566063dSJacob Faibussowitsch   PetscCall(GmshNodesCreate(num, &nodes));
5340598e1a2SLisandro Dalcin   mesh->numNodes = num;
5350598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
5360598e1a2SLisandro Dalcin   for (n = 0; n < num; ++n) {
5370598e1a2SLisandro Dalcin     double *xyz = nodes->xyz + n*3;
5389566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM));
5399566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE));
5409566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
5419566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
5420598e1a2SLisandro Dalcin     nodes->id[n] = nid;
54381a1af93SMatthew G. Knepley     nodes->tag[n] = -1;
54411c56cb0SLisandro Dalcin   }
54511c56cb0SLisandro Dalcin   PetscFunctionReturn(0);
54611c56cb0SLisandro Dalcin }
54711c56cb0SLisandro Dalcin 
548de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
549de78e4feSLisandro Dalcin    file contents multiple times to figure out the true number of cells and facets
550de78e4feSLisandro Dalcin    in the given mesh. To make this more efficient we read the file contents only
551de78e4feSLisandro Dalcin    once and store them in memory, while determining the true number of cells. */
5520598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v22(GmshFile* gmsh, GmshMesh *mesh)
55311c56cb0SLisandro Dalcin {
554698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
555698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
556698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
557de78e4feSLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
5580598e1a2SLisandro Dalcin   int            i, c, p, num, ibuf[1+4+1000], snum;
5590598e1a2SLisandro Dalcin   int            cellType, numElem, numVerts, numNodes, numTags;
560df032b82SMatthew G. Knepley   GmshElement   *elements;
5610598e1a2SLisandro Dalcin   PetscInt      *nodeMap = gmsh->nodeMap;
562df032b82SMatthew G. Knepley 
563df032b82SMatthew G. Knepley   PetscFunctionBegin;
5649566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
565698a79b9SLisandro Dalcin   snum = sscanf(line, "%d", &num);
56608401ef6SPierre Jolivet   PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
5679566063dSJacob Faibussowitsch   PetscCall(GmshElementsCreate(num, &elements));
5680598e1a2SLisandro Dalcin   mesh->numElems = num;
5690598e1a2SLisandro Dalcin   mesh->elements = elements;
570698a79b9SLisandro Dalcin   for (c = 0; c < num;) {
5719566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM));
5729566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3));
5730598e1a2SLisandro Dalcin 
5740598e1a2SLisandro Dalcin     cellType = binary ? ibuf[0] : ibuf[1];
5750598e1a2SLisandro Dalcin     numElem  = binary ? ibuf[1] : 1;
576df032b82SMatthew G. Knepley     numTags  = ibuf[2];
5770598e1a2SLisandro Dalcin 
5789566063dSJacob Faibussowitsch     PetscCall(GmshCellTypeCheck(cellType));
5790598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
5800598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
5810598e1a2SLisandro Dalcin 
582df032b82SMatthew G. Knepley     for (i = 0; i < numElem; ++i, ++c) {
5830598e1a2SLisandro Dalcin       GmshElement *element = elements + c;
5840598e1a2SLisandro Dalcin       const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off;
5850598e1a2SLisandro Dalcin       const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1;
5869566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, ibuf+off, nint, NULL, PETSC_ENUM));
5879566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(ibuf+off, PETSC_ENUM, nint));
5880598e1a2SLisandro Dalcin       element->id  = id[0];
5890598e1a2SLisandro Dalcin       element->dim = GmshCellMap[cellType].dim;
5900598e1a2SLisandro Dalcin       element->cellType = cellType;
5910598e1a2SLisandro Dalcin       element->numVerts = numVerts;
5920598e1a2SLisandro Dalcin       element->numNodes = numNodes;
5930598e1a2SLisandro Dalcin       element->numTags  = PetscMin(numTags, 4);
5949566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes));
5950598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
5960598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = tags[p];
597df032b82SMatthew G. Knepley     }
598df032b82SMatthew G. Knepley   }
599df032b82SMatthew G. Knepley   PetscFunctionReturn(0);
600df032b82SMatthew G. Knepley }
6017d282ae0SMichael Lange 
602de78e4feSLisandro Dalcin /*
603de78e4feSLisandro Dalcin $Entities
604de78e4feSLisandro Dalcin   numPoints(unsigned long) numCurves(unsigned long)
605de78e4feSLisandro Dalcin     numSurfaces(unsigned long) numVolumes(unsigned long)
606de78e4feSLisandro Dalcin   // points
607de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
608de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
609de78e4feSLisandro Dalcin     numPhysicals(unsigned long) phyisicalTag[...](int)
610de78e4feSLisandro Dalcin   ...
611de78e4feSLisandro Dalcin   // curves
612de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
613de78e4feSLisandro Dalcin      boxMaxX(double) boxMaxY(double) boxMaxZ(double)
614de78e4feSLisandro Dalcin      numPhysicals(unsigned long) physicalTag[...](int)
615de78e4feSLisandro Dalcin      numBREPVert(unsigned long) tagBREPVert[...](int)
616de78e4feSLisandro Dalcin   ...
617de78e4feSLisandro Dalcin   // surfaces
618de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
619de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
620de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
621de78e4feSLisandro Dalcin     numBREPCurve(unsigned long) tagBREPCurve[...](int)
622de78e4feSLisandro Dalcin   ...
623de78e4feSLisandro Dalcin   // volumes
624de78e4feSLisandro Dalcin   tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
625de78e4feSLisandro Dalcin     boxMaxX(double) boxMaxY(double) boxMaxZ(double)
626de78e4feSLisandro Dalcin     numPhysicals(unsigned long) physicalTag[...](int)
627de78e4feSLisandro Dalcin     numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int)
628de78e4feSLisandro Dalcin   ...
629de78e4feSLisandro Dalcin $EndEntities
630de78e4feSLisandro Dalcin */
6310598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh)
632de78e4feSLisandro Dalcin {
633698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
634698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
635698a79b9SLisandro Dalcin   long           index, num, lbuf[4];
636730356f6SLisandro Dalcin   int            dim, eid, numTags, *ibuf, t;
637698a79b9SLisandro Dalcin   PetscInt       count[4], i;
638a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
639de78e4feSLisandro Dalcin 
640de78e4feSLisandro Dalcin   PetscFunctionBegin;
6419566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG));
6429566063dSJacob Faibussowitsch   if (byteSwap) PetscCall(PetscByteSwap(lbuf, PETSC_LONG, 4));
643698a79b9SLisandro Dalcin   for (i = 0; i < 4; ++i) count[i] = lbuf[i];
6449566063dSJacob Faibussowitsch   PetscCall(GmshEntitiesCreate(count, &mesh->entities));
645de78e4feSLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
646730356f6SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
6479566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM));
6489566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(&eid, PETSC_ENUM, 1));
6499566063dSJacob Faibussowitsch       PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity));
6509566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE));
6519566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6));
6529566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG));
6539566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1));
6549566063dSJacob Faibussowitsch       PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf));
6559566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM));
6569566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num));
657de78e4feSLisandro Dalcin       entity->numTags = numTags = (int) PetscMin(num, 4);
658de78e4feSLisandro Dalcin       for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t];
659de78e4feSLisandro Dalcin       if (dim == 0) continue;
6609566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG));
6619566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1));
6629566063dSJacob Faibussowitsch       PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf));
6639566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM));
6649566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num));
665de78e4feSLisandro Dalcin     }
666de78e4feSLisandro Dalcin   }
667de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
668de78e4feSLisandro Dalcin }
669de78e4feSLisandro Dalcin 
670de78e4feSLisandro Dalcin /*
671de78e4feSLisandro Dalcin $Nodes
672de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numNodes(unsigned long)
673de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long)
674de78e4feSLisandro Dalcin     tag(int) x(double) y(double) z(double)
675de78e4feSLisandro Dalcin     ...
676de78e4feSLisandro Dalcin   ...
677de78e4feSLisandro Dalcin $EndNodes
678de78e4feSLisandro Dalcin */
6790598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh)
680de78e4feSLisandro Dalcin {
681698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
682698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
6830598e1a2SLisandro Dalcin   long           block, node, n, numEntityBlocks, numTotalNodes, numNodes;
684de78e4feSLisandro Dalcin   int            info[3], nid;
6850598e1a2SLisandro Dalcin   GmshNodes      *nodes;
686de78e4feSLisandro Dalcin 
687de78e4feSLisandro Dalcin   PetscFunctionBegin;
6889566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG));
6899566063dSJacob Faibussowitsch   if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1));
6909566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG));
6919566063dSJacob Faibussowitsch   if (byteSwap) PetscCall(PetscByteSwap(&numTotalNodes, PETSC_LONG, 1));
6929566063dSJacob Faibussowitsch   PetscCall(GmshNodesCreate(numTotalNodes, &nodes));
6930598e1a2SLisandro Dalcin   mesh->numNodes = numTotalNodes;
6940598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
6950598e1a2SLisandro Dalcin   for (n = 0, block = 0; block < numEntityBlocks; ++block) {
6969566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM));
6979566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG));
6989566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(&numNodes, PETSC_LONG, 1));
699698a79b9SLisandro Dalcin     if (gmsh->binary) {
700277f51e8SBarry Smith       size_t nbytes = sizeof(int) + 3*sizeof(double);
701277f51e8SBarry Smith       char   *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */
7029566063dSJacob Faibussowitsch       PetscCall(GmshBufferGet(gmsh, numNodes, nbytes, &cbuf));
7039566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR));
7040598e1a2SLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++n) {
705de78e4feSLisandro Dalcin         char   *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int);
7060598e1a2SLisandro Dalcin         double *xyz = nodes->xyz + n*3;
7079566063dSJacob Faibussowitsch         if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cnid, PETSC_ENUM, 1));
7089566063dSJacob Faibussowitsch         if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cxyz, PETSC_DOUBLE, 3));
7099566063dSJacob Faibussowitsch         PetscCall(PetscMemcpy(&nid, cnid, sizeof(int)));
7109566063dSJacob Faibussowitsch         PetscCall(PetscMemcpy(xyz, cxyz, 3*sizeof(double)));
7119566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
7129566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
7130598e1a2SLisandro Dalcin         nodes->id[n] = nid;
71481a1af93SMatthew G. Knepley         nodes->tag[n] = -1;
715de78e4feSLisandro Dalcin       }
716de78e4feSLisandro Dalcin     } else {
7170598e1a2SLisandro Dalcin       for (node = 0; node < numNodes; ++node, ++n) {
7180598e1a2SLisandro Dalcin         double *xyz = nodes->xyz + n*3;
7199566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM));
7209566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE));
7219566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
7229566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
7230598e1a2SLisandro Dalcin         nodes->id[n] = nid;
72481a1af93SMatthew G. Knepley         nodes->tag[n] = -1;
725de78e4feSLisandro Dalcin       }
726de78e4feSLisandro Dalcin     }
727de78e4feSLisandro Dalcin   }
728de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
729de78e4feSLisandro Dalcin }
730de78e4feSLisandro Dalcin 
731de78e4feSLisandro Dalcin /*
732de78e4feSLisandro Dalcin $Elements
733de78e4feSLisandro Dalcin   numEntityBlocks(unsigned long) numElements(unsigned long)
734de78e4feSLisandro Dalcin   tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long)
735de78e4feSLisandro Dalcin     tag(int) numVert[...](int)
736de78e4feSLisandro Dalcin     ...
737de78e4feSLisandro Dalcin   ...
738de78e4feSLisandro Dalcin $EndElements
739de78e4feSLisandro Dalcin */
7400598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh)
741de78e4feSLisandro Dalcin {
742698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
743698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
744de78e4feSLisandro Dalcin   long           c, block, numEntityBlocks, numTotalElements, elem, numElements;
745de78e4feSLisandro Dalcin   int            p, info[3], *ibuf = NULL;
7460598e1a2SLisandro Dalcin   int            eid, dim, cellType, numVerts, numNodes, numTags;
747a5ba3d71SLisandro Dalcin   GmshEntity     *entity = NULL;
748de78e4feSLisandro Dalcin   GmshElement    *elements;
7490598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
750de78e4feSLisandro Dalcin 
751de78e4feSLisandro Dalcin   PetscFunctionBegin;
7529566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG));
7539566063dSJacob Faibussowitsch   if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1));
7549566063dSJacob Faibussowitsch   PetscCall(PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG));
7559566063dSJacob Faibussowitsch   if (byteSwap) PetscCall(PetscByteSwap(&numTotalElements, PETSC_LONG, 1));
7569566063dSJacob Faibussowitsch   PetscCall(GmshElementsCreate(numTotalElements, &elements));
7570598e1a2SLisandro Dalcin   mesh->numElems = numTotalElements;
7580598e1a2SLisandro Dalcin   mesh->elements = elements;
759de78e4feSLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
7609566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM));
7619566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(info, PETSC_ENUM, 3));
762de78e4feSLisandro Dalcin     eid = info[0]; dim = info[1]; cellType = info[2];
7639566063dSJacob Faibussowitsch     PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
7649566063dSJacob Faibussowitsch     PetscCall(GmshCellTypeCheck(cellType));
7650598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
7660598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
767730356f6SLisandro Dalcin     numTags  = entity->numTags;
7689566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG));
7699566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(&numElements, PETSC_LONG, 1));
7709566063dSJacob Faibussowitsch     PetscCall(GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf));
7719566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM));
7729566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements));
773de78e4feSLisandro Dalcin     for (elem = 0; elem < numElements; ++elem, ++c) {
774de78e4feSLisandro Dalcin       GmshElement *element = elements + c;
7750598e1a2SLisandro Dalcin       const int *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
7760598e1a2SLisandro Dalcin       element->id  = id[0];
777de78e4feSLisandro Dalcin       element->dim = dim;
778de78e4feSLisandro Dalcin       element->cellType = cellType;
7790598e1a2SLisandro Dalcin       element->numVerts = numVerts;
780de78e4feSLisandro Dalcin       element->numNodes = numNodes;
781de78e4feSLisandro Dalcin       element->numTags  = numTags;
7829566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes));
7830598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
7840598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = entity->tags[p];
785de78e4feSLisandro Dalcin     }
786de78e4feSLisandro Dalcin   }
787698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
788698a79b9SLisandro Dalcin }
789698a79b9SLisandro Dalcin 
7900598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[])
791698a79b9SLisandro Dalcin {
792698a79b9SLisandro Dalcin   PetscViewer    viewer = gmsh->viewer;
793698a79b9SLisandro Dalcin   int            fileFormat = gmsh->fileFormat;
794698a79b9SLisandro Dalcin   PetscBool      binary = gmsh->binary;
795698a79b9SLisandro Dalcin   PetscBool      byteSwap = gmsh->byteSwap;
796698a79b9SLisandro Dalcin   int            numPeriodic, snum, i;
797698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
7980598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
799698a79b9SLisandro Dalcin 
800698a79b9SLisandro Dalcin   PetscFunctionBegin;
801698a79b9SLisandro Dalcin   if (fileFormat == 22 || !binary) {
8029566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
803698a79b9SLisandro Dalcin     snum = sscanf(line, "%d", &numPeriodic);
80408401ef6SPierre Jolivet     PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
805698a79b9SLisandro Dalcin   } else {
8069566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM));
8079566063dSJacob Faibussowitsch     if (byteSwap) PetscCall(PetscByteSwap(&numPeriodic, PETSC_ENUM, 1));
808698a79b9SLisandro Dalcin   }
809698a79b9SLisandro Dalcin   for (i = 0; i < numPeriodic; i++) {
8109dddd249SSatish Balay     int    ibuf[3], correspondingDim = -1, correspondingTag = -1, primaryTag = -1, correspondingNode, primaryNode;
811698a79b9SLisandro Dalcin     long   j, nNodes;
812698a79b9SLisandro Dalcin     double affine[16];
813698a79b9SLisandro Dalcin 
814698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
8159566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING));
8169dddd249SSatish Balay       snum = sscanf(line, "%d %d %d", &correspondingDim, &correspondingTag, &primaryTag);
81708401ef6SPierre Jolivet       PetscCheck(snum == 3,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
818698a79b9SLisandro Dalcin     } else {
8199566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM));
8209566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3));
8219dddd249SSatish Balay       correspondingDim = ibuf[0]; correspondingTag = ibuf[1]; primaryTag = ibuf[2];
822698a79b9SLisandro Dalcin     }
8239dddd249SSatish Balay     (void)correspondingDim; (void)correspondingTag; (void)primaryTag; /* unused */
824698a79b9SLisandro Dalcin 
825698a79b9SLisandro Dalcin     if (fileFormat == 22 || !binary) {
8269566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
827698a79b9SLisandro Dalcin       snum = sscanf(line, "%ld", &nNodes);
8284c500f23SPierre Jolivet       if (snum != 1) { /* discard transformation and try again */
8299566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING));
8309566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
831698a79b9SLisandro Dalcin         snum = sscanf(line, "%ld", &nNodes);
83208401ef6SPierre Jolivet         PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
833698a79b9SLisandro Dalcin       }
834698a79b9SLisandro Dalcin     } else {
8359566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG));
8369566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1));
8374c500f23SPierre Jolivet       if (nNodes == -1) { /* discard transformation and try again */
8389566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE));
8399566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG));
8409566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1));
841698a79b9SLisandro Dalcin       }
842698a79b9SLisandro Dalcin     }
843698a79b9SLisandro Dalcin 
844698a79b9SLisandro Dalcin     for (j = 0; j < nNodes; j++) {
845698a79b9SLisandro Dalcin       if (fileFormat == 22 || !binary) {
8469566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING));
8479dddd249SSatish Balay         snum = sscanf(line, "%d %d", &correspondingNode, &primaryNode);
84808401ef6SPierre Jolivet         PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
849698a79b9SLisandro Dalcin       } else {
8509566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM));
8519566063dSJacob Faibussowitsch         if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 2));
8529dddd249SSatish Balay         correspondingNode = ibuf[0]; primaryNode = ibuf[1];
853698a79b9SLisandro Dalcin       }
8549dddd249SSatish Balay       correspondingNode  = (int) nodeMap[correspondingNode];
8559dddd249SSatish Balay       primaryNode = (int) nodeMap[primaryNode];
8569dddd249SSatish Balay       periodicMap[correspondingNode] = primaryNode;
857698a79b9SLisandro Dalcin     }
858698a79b9SLisandro Dalcin   }
859698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
860698a79b9SLisandro Dalcin }
861698a79b9SLisandro Dalcin 
8620598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
863698a79b9SLisandro Dalcin $Entities
864698a79b9SLisandro Dalcin   numPoints(size_t) numCurves(size_t)
865698a79b9SLisandro Dalcin     numSurfaces(size_t) numVolumes(size_t)
866698a79b9SLisandro Dalcin   pointTag(int) X(double) Y(double) Z(double)
867698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
868698a79b9SLisandro Dalcin   ...
869698a79b9SLisandro Dalcin   curveTag(int) minX(double) minY(double) minZ(double)
870698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
871698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
872698a79b9SLisandro Dalcin     numBoundingPoints(size_t) pointTag(int) ...
873698a79b9SLisandro Dalcin   ...
874698a79b9SLisandro Dalcin   surfaceTag(int) minX(double) minY(double) minZ(double)
875698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
876698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
877698a79b9SLisandro Dalcin     numBoundingCurves(size_t) curveTag(int) ...
878698a79b9SLisandro Dalcin   ...
879698a79b9SLisandro Dalcin   volumeTag(int) minX(double) minY(double) minZ(double)
880698a79b9SLisandro Dalcin     maxX(double) maxY(double) maxZ(double)
881698a79b9SLisandro Dalcin     numPhysicalTags(size_t) physicalTag(int) ...
882698a79b9SLisandro Dalcin     numBoundngSurfaces(size_t) surfaceTag(int) ...
883698a79b9SLisandro Dalcin   ...
884698a79b9SLisandro Dalcin $EndEntities
885698a79b9SLisandro Dalcin */
8860598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh)
887698a79b9SLisandro Dalcin {
888698a79b9SLisandro Dalcin   PetscInt       count[4], index, numTags, i;
889698a79b9SLisandro Dalcin   int            dim, eid, *tags = NULL;
890698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
891698a79b9SLisandro Dalcin 
892698a79b9SLisandro Dalcin   PetscFunctionBegin;
8939566063dSJacob Faibussowitsch   PetscCall(GmshReadSize(gmsh, count, 4));
8949566063dSJacob Faibussowitsch   PetscCall(GmshEntitiesCreate(count, &mesh->entities));
895698a79b9SLisandro Dalcin   for (dim = 0; dim < 4; ++dim) {
896698a79b9SLisandro Dalcin     for (index = 0; index < count[dim]; ++index) {
8979566063dSJacob Faibussowitsch       PetscCall(GmshReadInt(gmsh, &eid, 1));
8989566063dSJacob Faibussowitsch       PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity));
8999566063dSJacob Faibussowitsch       PetscCall(GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6));
9009566063dSJacob Faibussowitsch       PetscCall(GmshReadSize(gmsh, &numTags, 1));
9019566063dSJacob Faibussowitsch       PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags));
9029566063dSJacob Faibussowitsch       PetscCall(GmshReadInt(gmsh, tags, numTags));
903698a79b9SLisandro Dalcin       entity->numTags = PetscMin(numTags, 4);
904698a79b9SLisandro Dalcin       for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i];
905698a79b9SLisandro Dalcin       if (dim == 0) continue;
9069566063dSJacob Faibussowitsch       PetscCall(GmshReadSize(gmsh, &numTags, 1));
9079566063dSJacob Faibussowitsch       PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags));
9089566063dSJacob Faibussowitsch       PetscCall(GmshReadInt(gmsh, tags, numTags));
90981a1af93SMatthew G. Knepley       /* Currently, we do not save the ids for the bounding entities */
910698a79b9SLisandro Dalcin     }
911698a79b9SLisandro Dalcin   }
912698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
913698a79b9SLisandro Dalcin }
914698a79b9SLisandro Dalcin 
91533a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
916698a79b9SLisandro Dalcin $Nodes
917698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numNodes(size_t)
918698a79b9SLisandro Dalcin     minNodeTag(size_t) maxNodeTag(size_t)
919698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t)
920698a79b9SLisandro Dalcin     nodeTag(size_t)
921698a79b9SLisandro Dalcin     ...
922698a79b9SLisandro Dalcin     x(double) y(double) z(double)
923698a79b9SLisandro Dalcin        < u(double; if parametric and entityDim = 1 or entityDim = 2) >
924698a79b9SLisandro Dalcin        < v(double; if parametric and entityDim = 2) >
925698a79b9SLisandro Dalcin     ...
926698a79b9SLisandro Dalcin   ...
927698a79b9SLisandro Dalcin $EndNodes
928698a79b9SLisandro Dalcin */
9290598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh)
930698a79b9SLisandro Dalcin {
93181a1af93SMatthew G. Knepley   int            info[3], dim, eid, parametric;
93281a1af93SMatthew G. Knepley   PetscInt       sizes[4], numEntityBlocks, numTags, numNodes, numNodesBlock = 0, block, node, n;
93381a1af93SMatthew G. Knepley   GmshEntity     *entity = NULL;
9340598e1a2SLisandro Dalcin   GmshNodes      *nodes;
935698a79b9SLisandro Dalcin 
936698a79b9SLisandro Dalcin   PetscFunctionBegin;
9379566063dSJacob Faibussowitsch   PetscCall(GmshReadSize(gmsh, sizes, 4));
9380598e1a2SLisandro Dalcin   numEntityBlocks = sizes[0]; numNodes = sizes[1];
9399566063dSJacob Faibussowitsch   PetscCall(GmshNodesCreate(numNodes, &nodes));
9400598e1a2SLisandro Dalcin   mesh->numNodes = numNodes;
9410598e1a2SLisandro Dalcin   mesh->nodelist = nodes;
942698a79b9SLisandro Dalcin   for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) {
9439566063dSJacob Faibussowitsch     PetscCall(GmshReadInt(gmsh, info, 3));
94481a1af93SMatthew G. Knepley     dim = info[0]; eid = info[1]; parametric = info[2];
9459566063dSJacob Faibussowitsch     PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
94681a1af93SMatthew G. Knepley     numTags = PetscMin(1, entity->numTags);
94781a1af93SMatthew G. Knepley     if (entity->numTags > 1) PetscInfo(NULL, "Entity %d has more than %d physical tags, assigning only the first to nodes", eid, 1);
94881a1af93SMatthew G. Knepley     PetscCheck(!parametric, PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported");
9499566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, &numNodesBlock, 1));
9509566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, nodes->id+node, numNodesBlock));
9519566063dSJacob Faibussowitsch     PetscCall(GmshReadDouble(gmsh, nodes->xyz+node*3, numNodesBlock*3));
95281a1af93SMatthew G. Knepley     for (n = 0; n < numNodesBlock; ++n) nodes->tag[node+n] = numTags ? entity->tags[0] : -1;
953698a79b9SLisandro Dalcin   }
9540598e1a2SLisandro Dalcin   gmsh->nodeStart = sizes[2];
9550598e1a2SLisandro Dalcin   gmsh->nodeEnd   = sizes[3]+1;
956698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
957698a79b9SLisandro Dalcin }
958698a79b9SLisandro Dalcin 
95933a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
960698a79b9SLisandro Dalcin $Elements
961698a79b9SLisandro Dalcin   numEntityBlocks(size_t) numElements(size_t)
962698a79b9SLisandro Dalcin     minElementTag(size_t) maxElementTag(size_t)
963698a79b9SLisandro Dalcin   entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t)
964698a79b9SLisandro Dalcin     elementTag(size_t) nodeTag(size_t) ...
965698a79b9SLisandro Dalcin     ...
966698a79b9SLisandro Dalcin   ...
967698a79b9SLisandro Dalcin $EndElements
968698a79b9SLisandro Dalcin */
9690598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh)
970698a79b9SLisandro Dalcin {
9710598e1a2SLisandro Dalcin   int            info[3], eid, dim, cellType;
9720598e1a2SLisandro Dalcin   PetscInt       sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p;
973698a79b9SLisandro Dalcin   GmshEntity     *entity = NULL;
974698a79b9SLisandro Dalcin   GmshElement    *elements;
9750598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
976698a79b9SLisandro Dalcin 
977698a79b9SLisandro Dalcin   PetscFunctionBegin;
9789566063dSJacob Faibussowitsch   PetscCall(GmshReadSize(gmsh, sizes, 4));
979698a79b9SLisandro Dalcin   numEntityBlocks = sizes[0]; numElements = sizes[1];
9809566063dSJacob Faibussowitsch   PetscCall(GmshElementsCreate(numElements, &elements));
9810598e1a2SLisandro Dalcin   mesh->numElems = numElements;
9820598e1a2SLisandro Dalcin   mesh->elements = elements;
983698a79b9SLisandro Dalcin   for (c = 0, block = 0; block < numEntityBlocks; ++block) {
9849566063dSJacob Faibussowitsch     PetscCall(GmshReadInt(gmsh, info, 3));
985698a79b9SLisandro Dalcin     dim = info[0]; eid = info[1]; cellType = info[2];
9869566063dSJacob Faibussowitsch     PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
9879566063dSJacob Faibussowitsch     PetscCall(GmshCellTypeCheck(cellType));
9880598e1a2SLisandro Dalcin     numVerts = GmshCellMap[cellType].numVerts;
9890598e1a2SLisandro Dalcin     numNodes = GmshCellMap[cellType].numNodes;
99081a1af93SMatthew G. Knepley     numTags  = PetscMin(4, entity->numTags);
99181a1af93SMatthew G. Knepley     if (entity->numTags > 4) PetscInfo(NULL, "Entity %d has more then %d physical tags, assigning only the first to elements", eid, 4);
9929566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, &numBlockElements, 1));
9939566063dSJacob Faibussowitsch     PetscCall(GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf));
9949566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements));
995698a79b9SLisandro Dalcin     for (elem = 0; elem < numBlockElements; ++elem, ++c) {
996698a79b9SLisandro Dalcin       GmshElement *element = elements + c;
9970598e1a2SLisandro Dalcin       const PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1;
998698a79b9SLisandro Dalcin       element->id  = id[0];
999698a79b9SLisandro Dalcin       element->dim = dim;
1000698a79b9SLisandro Dalcin       element->cellType = cellType;
10010598e1a2SLisandro Dalcin       element->numVerts = numVerts;
1002698a79b9SLisandro Dalcin       element->numNodes = numNodes;
1003698a79b9SLisandro Dalcin       element->numTags  = numTags;
10049566063dSJacob Faibussowitsch       PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes));
10050598e1a2SLisandro Dalcin       for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
10060598e1a2SLisandro Dalcin       for (p = 0; p < element->numTags;  p++) element->tags[p]  = entity->tags[p];
1007698a79b9SLisandro Dalcin     }
1008698a79b9SLisandro Dalcin   }
1009698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
1010698a79b9SLisandro Dalcin }
1011698a79b9SLisandro Dalcin 
10120598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1013698a79b9SLisandro Dalcin $Periodic
1014698a79b9SLisandro Dalcin   numPeriodicLinks(size_t)
10159dddd249SSatish Balay   entityDim(int) entityTag(int) entityTagPrimary(int)
1016698a79b9SLisandro Dalcin   numAffine(size_t) value(double) ...
1017698a79b9SLisandro Dalcin   numCorrespondingNodes(size_t)
10189dddd249SSatish Balay     nodeTag(size_t) nodeTagPrimary(size_t)
1019698a79b9SLisandro Dalcin     ...
1020698a79b9SLisandro Dalcin   ...
1021698a79b9SLisandro Dalcin $EndPeriodic
1022698a79b9SLisandro Dalcin */
10230598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[])
1024698a79b9SLisandro Dalcin {
1025698a79b9SLisandro Dalcin   int            info[3];
1026698a79b9SLisandro Dalcin   double         dbuf[16];
10270598e1a2SLisandro Dalcin   PetscInt       numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node;
10280598e1a2SLisandro Dalcin   PetscInt       *nodeMap = gmsh->nodeMap;
1029698a79b9SLisandro Dalcin 
1030698a79b9SLisandro Dalcin   PetscFunctionBegin;
10319566063dSJacob Faibussowitsch   PetscCall(GmshReadSize(gmsh, &numPeriodicLinks, 1));
1032698a79b9SLisandro Dalcin   for (link = 0; link < numPeriodicLinks; ++link) {
10339566063dSJacob Faibussowitsch     PetscCall(GmshReadInt(gmsh, info, 3));
10349566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, &numAffine, 1));
10359566063dSJacob Faibussowitsch     PetscCall(GmshReadDouble(gmsh, dbuf, numAffine));
10369566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, &numCorrespondingNodes, 1));
10379566063dSJacob Faibussowitsch     PetscCall(GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags));
10389566063dSJacob Faibussowitsch     PetscCall(GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2));
1039698a79b9SLisandro Dalcin     for (node = 0; node < numCorrespondingNodes; ++node) {
10409dddd249SSatish Balay       PetscInt correspondingNode = nodeMap[nodeTags[node*2+0]];
10419dddd249SSatish Balay       PetscInt primaryNode = nodeMap[nodeTags[node*2+1]];
10429dddd249SSatish Balay       periodicMap[correspondingNode] = primaryNode;
1043698a79b9SLisandro Dalcin     }
1044698a79b9SLisandro Dalcin   }
1045698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
1046698a79b9SLisandro Dalcin }
1047698a79b9SLisandro Dalcin 
10480598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1049d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2
1050d6689ca9SLisandro Dalcin   version(ASCII double; currently 4.1)
1051d6689ca9SLisandro Dalcin   file-type(ASCII int; 0 for ASCII mode, 1 for binary mode)
1052d6689ca9SLisandro Dalcin   data-size(ASCII int; sizeof(size_t))
1053d6689ca9SLisandro Dalcin   < int with value one; only in binary mode, to detect endianness >
1054d6689ca9SLisandro Dalcin $EndMeshFormat
1055d6689ca9SLisandro Dalcin */
10560598e1a2SLisandro Dalcin static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh)
1057698a79b9SLisandro Dalcin {
1058698a79b9SLisandro Dalcin   char           line[PETSC_MAX_PATH_LEN];
1059698a79b9SLisandro Dalcin   int            snum, fileType, fileFormat, dataSize, checkEndian;
1060698a79b9SLisandro Dalcin   float          version;
1061698a79b9SLisandro Dalcin 
1062698a79b9SLisandro Dalcin   PetscFunctionBegin;
10639566063dSJacob Faibussowitsch   PetscCall(GmshReadString(gmsh, line, 3));
1064698a79b9SLisandro Dalcin   snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
106508401ef6SPierre Jolivet   PetscCheck(snum == 3,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
106608401ef6SPierre Jolivet   PetscCheck(version >= 2.2,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
10672c71b3e2SJacob Faibussowitsch   PetscCheckFalse((int)version == 3,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
106808401ef6SPierre Jolivet   PetscCheck(version <= 4.1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
106908401ef6SPierre Jolivet   PetscCheck(!gmsh->binary || fileType,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII");
107008401ef6SPierre Jolivet   PetscCheck(gmsh->binary || !fileType,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary");
1071af7a0b9fSSatish Balay   fileFormat = (int)roundf(version*10);
10722c71b3e2SJacob Faibussowitsch   PetscCheckFalse(fileFormat <= 40 && dataSize != sizeof(double),PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize);
10732c71b3e2SJacob Faibussowitsch   PetscCheckFalse(fileFormat >= 41 && dataSize != sizeof(int) && dataSize != sizeof(PetscInt64),PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize);
1074698a79b9SLisandro Dalcin   gmsh->fileFormat = fileFormat;
1075698a79b9SLisandro Dalcin   gmsh->dataSize = dataSize;
1076698a79b9SLisandro Dalcin   gmsh->byteSwap = PETSC_FALSE;
1077698a79b9SLisandro Dalcin   if (gmsh->binary) {
10789566063dSJacob Faibussowitsch     PetscCall(GmshReadInt(gmsh, &checkEndian, 1));
1079698a79b9SLisandro Dalcin     if (checkEndian != 1) {
10809566063dSJacob Faibussowitsch       PetscCall(PetscByteSwap(&checkEndian, PETSC_ENUM, 1));
108108401ef6SPierre Jolivet       PetscCheck(checkEndian == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line);
1082698a79b9SLisandro Dalcin       gmsh->byteSwap = PETSC_TRUE;
1083698a79b9SLisandro Dalcin     }
1084698a79b9SLisandro Dalcin   }
1085698a79b9SLisandro Dalcin   PetscFunctionReturn(0);
1086698a79b9SLisandro Dalcin }
1087698a79b9SLisandro Dalcin 
10881f643af3SLisandro Dalcin /*
10891f643af3SLisandro Dalcin PhysicalNames
10901f643af3SLisandro Dalcin   numPhysicalNames(ASCII int)
10911f643af3SLisandro Dalcin   dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max)
10921f643af3SLisandro Dalcin   ...
10931f643af3SLisandro Dalcin $EndPhysicalNames
10941f643af3SLisandro Dalcin */
1095a45dabc8SMatthew G. Knepley static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh, GmshMesh *mesh)
1096698a79b9SLisandro Dalcin {
1097*5f5cd6d5SMatthew G. Knepley   char           line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q, *r;
1098a45dabc8SMatthew G. Knepley   int            snum, region, dim, tag;
1099698a79b9SLisandro Dalcin 
1100698a79b9SLisandro Dalcin   PetscFunctionBegin;
11019566063dSJacob Faibussowitsch   PetscCall(GmshReadString(gmsh, line, 1));
1102a45dabc8SMatthew G. Knepley   snum = sscanf(line, "%d", &region);
1103a45dabc8SMatthew G. Knepley   mesh->numRegions = region;
110408401ef6SPierre Jolivet   PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
11059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(mesh->numRegions, &mesh->regionTags, mesh->numRegions, &mesh->regionNames));
1106a45dabc8SMatthew G. Knepley   for (region = 0; region < mesh->numRegions; ++region) {
11079566063dSJacob Faibussowitsch     PetscCall(GmshReadString(gmsh, line, 2));
11081f643af3SLisandro Dalcin     snum = sscanf(line, "%d %d", &dim, &tag);
110908401ef6SPierre Jolivet     PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
11109566063dSJacob Faibussowitsch     PetscCall(GmshReadString(gmsh, line, -(PetscInt)sizeof(line)));
11119566063dSJacob Faibussowitsch     PetscCall(PetscStrchr(line, '"', &p));
111228b400f6SJacob Faibussowitsch     PetscCheck(p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
11139566063dSJacob Faibussowitsch     PetscCall(PetscStrrchr(line, '"', &q));
111408401ef6SPierre Jolivet     PetscCheck(q != p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1115*5f5cd6d5SMatthew G. Knepley     PetscCall(PetscStrrchr(line, ':', &r));
1116*5f5cd6d5SMatthew G. Knepley     if (p != r) q = r;
11179566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(name, p+1, (size_t)(q-p-1)));
1118a45dabc8SMatthew G. Knepley     mesh->regionTags[region] = tag;
11199566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(name, &mesh->regionNames[region]));
1120698a79b9SLisandro Dalcin   }
1121de78e4feSLisandro Dalcin   PetscFunctionReturn(0);
1122de78e4feSLisandro Dalcin }
1123de78e4feSLisandro Dalcin 
11240598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh)
112596ca5757SLisandro Dalcin {
11260598e1a2SLisandro Dalcin   PetscFunctionBegin;
11270598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11289566063dSJacob Faibussowitsch   case 41: PetscCall(GmshReadEntities_v41(gmsh, mesh)); break;
11299566063dSJacob Faibussowitsch   default: PetscCall(GmshReadEntities_v40(gmsh, mesh)); break;
113096ca5757SLisandro Dalcin   }
11310598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
11320598e1a2SLisandro Dalcin }
11330598e1a2SLisandro Dalcin 
11340598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh)
11350598e1a2SLisandro Dalcin {
11360598e1a2SLisandro Dalcin   PetscFunctionBegin;
11370598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11389566063dSJacob Faibussowitsch   case 41: PetscCall(GmshReadNodes_v41(gmsh, mesh)); break;
11399566063dSJacob Faibussowitsch   case 40: PetscCall(GmshReadNodes_v40(gmsh, mesh)); break;
11409566063dSJacob Faibussowitsch   default: PetscCall(GmshReadNodes_v22(gmsh, mesh)); break;
11410598e1a2SLisandro Dalcin   }
11420598e1a2SLisandro Dalcin 
11430598e1a2SLisandro Dalcin   { /* Gmsh v2.2/v4.0 does not provide min/max node tags */
11440598e1a2SLisandro Dalcin     if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) {
11450598e1a2SLisandro Dalcin       PetscInt  tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n;
11460598e1a2SLisandro Dalcin       GmshNodes *nodes = mesh->nodelist;
11470598e1a2SLisandro Dalcin       for (n = 0; n < mesh->numNodes; ++n) {
11480598e1a2SLisandro Dalcin         const PetscInt tag = nodes->id[n];
11490598e1a2SLisandro Dalcin         tagMin = PetscMin(tag, tagMin);
11500598e1a2SLisandro Dalcin         tagMax = PetscMax(tag, tagMax);
11510598e1a2SLisandro Dalcin       }
11520598e1a2SLisandro Dalcin       gmsh->nodeStart = tagMin;
11530598e1a2SLisandro Dalcin       gmsh->nodeEnd   = tagMax+1;
11540598e1a2SLisandro Dalcin     }
11550598e1a2SLisandro Dalcin   }
11560598e1a2SLisandro Dalcin 
11570598e1a2SLisandro Dalcin   { /* Support for sparse node tags */
11580598e1a2SLisandro Dalcin     PetscInt  n, t;
11590598e1a2SLisandro Dalcin     GmshNodes *nodes = mesh->nodelist;
11609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf));
11610598e1a2SLisandro Dalcin     for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT;
11620598e1a2SLisandro Dalcin     gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart;
11630598e1a2SLisandro Dalcin     for (n = 0; n < mesh->numNodes; ++n) {
11640598e1a2SLisandro Dalcin       const PetscInt tag = nodes->id[n];
116508401ef6SPierre Jolivet       PetscCheck(gmsh->nodeMap[tag] < 0,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %D", tag);
11660598e1a2SLisandro Dalcin       gmsh->nodeMap[tag] = n;
11670598e1a2SLisandro Dalcin     }
11680598e1a2SLisandro Dalcin   }
11690598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
11700598e1a2SLisandro Dalcin }
11710598e1a2SLisandro Dalcin 
11720598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh)
11730598e1a2SLisandro Dalcin {
11740598e1a2SLisandro Dalcin   PetscFunctionBegin;
11750598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
11769566063dSJacob Faibussowitsch   case 41: PetscCall(GmshReadElements_v41(gmsh, mesh)); break;
11779566063dSJacob Faibussowitsch   case 40: PetscCall(GmshReadElements_v40(gmsh, mesh)); break;
11789566063dSJacob Faibussowitsch   default: PetscCall(GmshReadElements_v22(gmsh, mesh)); break;
11790598e1a2SLisandro Dalcin   }
11800598e1a2SLisandro Dalcin 
11810598e1a2SLisandro Dalcin   { /* Reorder elements by codimension and polytope type */
11820598e1a2SLisandro Dalcin     PetscInt    ne = mesh->numElems;
11830598e1a2SLisandro Dalcin     GmshElement *elements = mesh->elements;
1184066ea43fSLisandro Dalcin     PetscInt    keymap[GMSH_NUM_POLYTOPES], nk = 0;
1185066ea43fSLisandro Dalcin     PetscInt    offset[GMSH_NUM_POLYTOPES+1], e, k;
11860598e1a2SLisandro Dalcin 
1187066ea43fSLisandro Dalcin     for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT;
11889566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(offset,sizeof(offset)));
11890598e1a2SLisandro Dalcin 
1190066ea43fSLisandro Dalcin     keymap[GMSH_TET] = nk++;
1191066ea43fSLisandro Dalcin     keymap[GMSH_HEX] = nk++;
1192066ea43fSLisandro Dalcin     keymap[GMSH_PRI] = nk++;
1193066ea43fSLisandro Dalcin     keymap[GMSH_PYR] = nk++;
1194066ea43fSLisandro Dalcin     keymap[GMSH_TRI] = nk++;
1195066ea43fSLisandro Dalcin     keymap[GMSH_QUA] = nk++;
1196066ea43fSLisandro Dalcin     keymap[GMSH_SEG] = nk++;
1197066ea43fSLisandro Dalcin     keymap[GMSH_VTX] = nk++;
11980598e1a2SLisandro Dalcin 
11999566063dSJacob Faibussowitsch     PetscCall(GmshElementsCreate(mesh->numElems, &mesh->elements));
12000598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope]
12010598e1a2SLisandro Dalcin     for (e = 0; e < ne; ++e) offset[1+key(e)]++;
12020598e1a2SLisandro Dalcin     for (k = 1; k < nk; ++k) offset[k] += offset[k-1];
12030598e1a2SLisandro Dalcin     for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e];
12040598e1a2SLisandro Dalcin #undef key
12059566063dSJacob Faibussowitsch     PetscCall(GmshElementsDestroy(&elements));
1206066ea43fSLisandro Dalcin   }
12070598e1a2SLisandro Dalcin 
1208066ea43fSLisandro Dalcin   { /* Mesh dimension and order */
1209066ea43fSLisandro Dalcin     GmshElement *elem = mesh->numElems ? mesh->elements : NULL;
1210066ea43fSLisandro Dalcin     mesh->dim   = elem ? GmshCellMap[elem->cellType].dim   : 0;
1211066ea43fSLisandro Dalcin     mesh->order = elem ? GmshCellMap[elem->cellType].order : 0;
12120598e1a2SLisandro Dalcin   }
12130598e1a2SLisandro Dalcin 
12140598e1a2SLisandro Dalcin   {
12150598e1a2SLisandro Dalcin     PetscBT  vtx;
12160598e1a2SLisandro Dalcin     PetscInt dim = mesh->dim, e, n, v;
12170598e1a2SLisandro Dalcin 
12189566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(mesh->numNodes, &vtx));
12190598e1a2SLisandro Dalcin 
12200598e1a2SLisandro Dalcin     /* Compute number of cells and set of vertices */
12210598e1a2SLisandro Dalcin     mesh->numCells = 0;
12220598e1a2SLisandro Dalcin     for (e = 0; e < mesh->numElems; ++e) {
12230598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + e;
12240598e1a2SLisandro Dalcin       if (elem->dim == dim && dim > 0) mesh->numCells++;
12250598e1a2SLisandro Dalcin       for (v = 0; v < elem->numVerts; v++) {
12269566063dSJacob Faibussowitsch         PetscCall(PetscBTSet(vtx, elem->nodes[v]));
12270598e1a2SLisandro Dalcin       }
12280598e1a2SLisandro Dalcin     }
12290598e1a2SLisandro Dalcin 
12300598e1a2SLisandro Dalcin     /* Compute numbering for vertices */
12310598e1a2SLisandro Dalcin     mesh->numVerts = 0;
12329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(mesh->numNodes, &mesh->vertexMap));
12330598e1a2SLisandro Dalcin     for (n = 0; n < mesh->numNodes; ++n)
12340598e1a2SLisandro Dalcin       mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT;
12350598e1a2SLisandro Dalcin 
12369566063dSJacob Faibussowitsch     PetscCall(PetscBTDestroy(&vtx));
12370598e1a2SLisandro Dalcin   }
12380598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
12390598e1a2SLisandro Dalcin }
12400598e1a2SLisandro Dalcin 
12410598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh)
12420598e1a2SLisandro Dalcin {
12430598e1a2SLisandro Dalcin   PetscInt       n;
12440598e1a2SLisandro Dalcin 
12450598e1a2SLisandro Dalcin   PetscFunctionBegin;
12469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(mesh->numNodes, &mesh->periodMap));
12470598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n;
12480598e1a2SLisandro Dalcin   switch (gmsh->fileFormat) {
12499566063dSJacob Faibussowitsch   case 41: PetscCall(GmshReadPeriodic_v41(gmsh, mesh->periodMap)); break;
12509566063dSJacob Faibussowitsch   default: PetscCall(GmshReadPeriodic_v40(gmsh, mesh->periodMap)); break;
12510598e1a2SLisandro Dalcin   }
12520598e1a2SLisandro Dalcin 
12539dddd249SSatish Balay   /* Find canonical primary nodes */
12540598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12550598e1a2SLisandro Dalcin     while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]])
12560598e1a2SLisandro Dalcin       mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]];
12570598e1a2SLisandro Dalcin 
12589dddd249SSatish Balay   /* Renumber vertices (filter out correspondings) */
12590598e1a2SLisandro Dalcin   mesh->numVerts = 0;
12600598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12610598e1a2SLisandro Dalcin     if (mesh->vertexMap[n] >= 0)   /* is vertex */
12629dddd249SSatish Balay       if (mesh->periodMap[n] == n) /* is primary */
12630598e1a2SLisandro Dalcin         mesh->vertexMap[n] = mesh->numVerts++;
12640598e1a2SLisandro Dalcin   for (n = 0; n < mesh->numNodes; ++n)
12650598e1a2SLisandro Dalcin     if (mesh->vertexMap[n] >= 0)   /* is vertex */
12669dddd249SSatish Balay       if (mesh->periodMap[n] != n) /* is corresponding  */
12670598e1a2SLisandro Dalcin         mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]];
12680598e1a2SLisandro Dalcin   PetscFunctionReturn(0);
12690598e1a2SLisandro Dalcin }
12700598e1a2SLisandro Dalcin 
1271066ea43fSLisandro Dalcin #define DM_POLYTOPE_VERTEX  DM_POLYTOPE_POINT
1272066ea43fSLisandro Dalcin #define DM_POLYTOPE_PYRAMID DM_POLYTOPE_UNKNOWN
1273066ea43fSLisandro Dalcin static const DMPolytopeType DMPolytopeMap[] = {
1274066ea43fSLisandro Dalcin   /* GMSH_VTX */ DM_POLYTOPE_VERTEX,
1275066ea43fSLisandro Dalcin   /* GMSH_SEG */ DM_POLYTOPE_SEGMENT,
1276066ea43fSLisandro Dalcin   /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE,
1277066ea43fSLisandro Dalcin   /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL,
1278066ea43fSLisandro Dalcin   /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON,
1279066ea43fSLisandro Dalcin   /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON,
1280066ea43fSLisandro Dalcin   /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM,
1281066ea43fSLisandro Dalcin   /* GMSH_PYR */ DM_POLYTOPE_PYRAMID,
1282066ea43fSLisandro Dalcin   DM_POLYTOPE_UNKNOWN
1283066ea43fSLisandro Dalcin };
12840598e1a2SLisandro Dalcin 
12859fbee547SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType)
12860598e1a2SLisandro Dalcin {
1287066ea43fSLisandro Dalcin   return DMPolytopeMap[GmshCellMap[cellType].polytope];
1288066ea43fSLisandro Dalcin }
1289066ea43fSLisandro Dalcin 
1290066ea43fSLisandro Dalcin static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem)
1291066ea43fSLisandro Dalcin {
1292066ea43fSLisandro Dalcin   DM              K;
1293066ea43fSLisandro Dalcin   PetscSpace      P;
1294066ea43fSLisandro Dalcin   PetscDualSpace  Q;
1295066ea43fSLisandro Dalcin   PetscQuadrature q, fq;
1296066ea43fSLisandro Dalcin   PetscBool       isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1297066ea43fSLisandro Dalcin   PetscBool       endpoint = PETSC_TRUE;
1298066ea43fSLisandro Dalcin   char            name[32];
1299066ea43fSLisandro Dalcin 
1300066ea43fSLisandro Dalcin   PetscFunctionBegin;
1301066ea43fSLisandro Dalcin   /* Create space */
13029566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(comm, &P));
13039566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
13049566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialSetTensor(P, isTensor));
13059566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(P, Nc));
13069566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(P, dim));
13079566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetDegree(P, k, PETSC_DETERMINE));
1308066ea43fSLisandro Dalcin   if (prefix) {
13099566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix));
13109566063dSJacob Faibussowitsch     PetscCall(PetscSpaceSetFromOptions(P));
13119566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, NULL));
13129566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetDegree(P, &k, NULL));
1313066ea43fSLisandro Dalcin   }
13149566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetUp(P));
1315066ea43fSLisandro Dalcin   /* Create dual space */
13169566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceCreate(comm, &Q));
13179566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
13189566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor));
13199566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity));
13209566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0));
13219566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
13229566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetOrder(Q, k));
13239566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
13249566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetDM(Q, K));
13259566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&K));
1326066ea43fSLisandro Dalcin   if (prefix) {
13279566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix));
13289566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceSetFromOptions(Q));
13299566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, NULL));
1330066ea43fSLisandro Dalcin   }
13319566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceSetUp(Q));
1332066ea43fSLisandro Dalcin   /* Create quadrature */
1333066ea43fSLisandro Dalcin   if (isSimplex) {
13349566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim,   1, k+1, -1, +1, &q));
13359566063dSJacob Faibussowitsch     PetscCall(PetscDTStroudConicalQuadrature(dim-1, 1, k+1, -1, +1, &fq));
1336066ea43fSLisandro Dalcin   } else {
13379566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim,   1, k+1, -1, +1, &q));
13389566063dSJacob Faibussowitsch     PetscCall(PetscDTGaussTensorQuadrature(dim-1, 1, k+1, -1, +1, &fq));
1339066ea43fSLisandro Dalcin   }
1340066ea43fSLisandro Dalcin   /* Create finite element */
13419566063dSJacob Faibussowitsch   PetscCall(PetscFECreate(comm, fem));
13429566063dSJacob Faibussowitsch   PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
13439566063dSJacob Faibussowitsch   PetscCall(PetscFESetNumComponents(*fem, Nc));
13449566063dSJacob Faibussowitsch   PetscCall(PetscFESetBasisSpace(*fem, P));
13459566063dSJacob Faibussowitsch   PetscCall(PetscFESetDualSpace(*fem, Q));
13469566063dSJacob Faibussowitsch   PetscCall(PetscFESetQuadrature(*fem, q));
13479566063dSJacob Faibussowitsch   PetscCall(PetscFESetFaceQuadrature(*fem, fq));
1348066ea43fSLisandro Dalcin   if (prefix) {
13499566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix));
13509566063dSJacob Faibussowitsch     PetscCall(PetscFESetFromOptions(*fem));
13519566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, NULL));
1352066ea43fSLisandro Dalcin   }
13539566063dSJacob Faibussowitsch   PetscCall(PetscFESetUp(*fem));
1354066ea43fSLisandro Dalcin   /* Cleanup */
13559566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&P));
13569566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&Q));
13579566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&q));
13589566063dSJacob Faibussowitsch   PetscCall(PetscQuadratureDestroy(&fq));
1359066ea43fSLisandro Dalcin   /* Set finite element name */
13609566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(name, sizeof(name), "%s%D", isSimplex? "P" : "Q", k));
13619566063dSJacob Faibussowitsch   PetscCall(PetscFESetName(*fem, name));
1362066ea43fSLisandro Dalcin   PetscFunctionReturn(0);
136396ca5757SLisandro Dalcin }
136496ca5757SLisandro Dalcin 
1365d6689ca9SLisandro Dalcin /*@C
1366d6689ca9SLisandro Dalcin   DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file
1367d6689ca9SLisandro Dalcin 
1368d6689ca9SLisandro Dalcin + comm        - The MPI communicator
1369d6689ca9SLisandro Dalcin . filename    - Name of the Gmsh file
1370d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh
1371d6689ca9SLisandro Dalcin 
1372d6689ca9SLisandro Dalcin   Output Parameter:
1373d6689ca9SLisandro Dalcin . dm  - The DM object representing the mesh
1374d6689ca9SLisandro Dalcin 
1375d6689ca9SLisandro Dalcin   Level: beginner
1376d6689ca9SLisandro Dalcin 
1377d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate()
1378d6689ca9SLisandro Dalcin @*/
1379d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1380d6689ca9SLisandro Dalcin {
1381d6689ca9SLisandro Dalcin   PetscViewer     viewer;
1382d6689ca9SLisandro Dalcin   PetscMPIInt     rank;
1383d6689ca9SLisandro Dalcin   int             fileType;
1384d6689ca9SLisandro Dalcin   PetscViewerType vtype;
1385d6689ca9SLisandro Dalcin 
1386d6689ca9SLisandro Dalcin   PetscFunctionBegin;
13879566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1388d6689ca9SLisandro Dalcin 
1389d6689ca9SLisandro Dalcin   /* Determine Gmsh file type (ASCII or binary) from file header */
1390dd400576SPatrick Sanan   if (rank == 0) {
13910598e1a2SLisandro Dalcin     GmshFile    gmsh[1];
1392d6689ca9SLisandro Dalcin     char        line[PETSC_MAX_PATH_LEN];
1393d6689ca9SLisandro Dalcin     int         snum;
1394d6689ca9SLisandro Dalcin     float       version;
1395d6689ca9SLisandro Dalcin 
13969566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(gmsh,1));
13979566063dSJacob Faibussowitsch     PetscCall(PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer));
13989566063dSJacob Faibussowitsch     PetscCall(PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII));
13999566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ));
14009566063dSJacob Faibussowitsch     PetscCall(PetscViewerFileSetName(gmsh->viewer, filename));
1401d6689ca9SLisandro Dalcin     /* Read only the first two lines of the Gmsh file */
14029566063dSJacob Faibussowitsch     PetscCall(GmshReadSection(gmsh, line));
14039566063dSJacob Faibussowitsch     PetscCall(GmshExpect(gmsh, "$MeshFormat", line));
14049566063dSJacob Faibussowitsch     PetscCall(GmshReadString(gmsh, line, 2));
1405d6689ca9SLisandro Dalcin     snum = sscanf(line, "%f %d", &version, &fileType);
140608401ef6SPierre Jolivet     PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
140708401ef6SPierre Jolivet     PetscCheck(version >= 2.2,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
14082c71b3e2SJacob Faibussowitsch     PetscCheckFalse((int)version == 3,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
140908401ef6SPierre Jolivet     PetscCheck(version <= 4.1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
14109566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&gmsh->viewer));
1411d6689ca9SLisandro Dalcin   }
14129566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&fileType, 1, MPI_INT, 0, comm));
1413d6689ca9SLisandro Dalcin   vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
1414d6689ca9SLisandro Dalcin 
1415d6689ca9SLisandro Dalcin   /* Create appropriate viewer and build plex */
14169566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, &viewer));
14179566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(viewer, vtype));
14189566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
14199566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(viewer, filename));
14209566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateGmsh(comm, viewer, interpolate, dm));
14219566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
1422d6689ca9SLisandro Dalcin   PetscFunctionReturn(0);
1423d6689ca9SLisandro Dalcin }
1424d6689ca9SLisandro Dalcin 
1425331830f3SMatthew G. Knepley /*@
14267d282ae0SMichael Lange   DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer
1427331830f3SMatthew G. Knepley 
1428d083f849SBarry Smith   Collective
1429331830f3SMatthew G. Knepley 
1430331830f3SMatthew G. Knepley   Input Parameters:
1431331830f3SMatthew G. Knepley + comm  - The MPI communicator
1432331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file
1433331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
1434331830f3SMatthew G. Knepley 
1435331830f3SMatthew G. Knepley   Output Parameter:
1436331830f3SMatthew G. Knepley . dm  - The DM object representing the mesh
1437331830f3SMatthew G. Knepley 
1438698a79b9SLisandro Dalcin   Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format
1439331830f3SMatthew G. Knepley 
1440331830f3SMatthew G. Knepley   Level: beginner
1441331830f3SMatthew G. Knepley 
1442331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate()
1443331830f3SMatthew G. Knepley @*/
1444331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
1445331830f3SMatthew G. Knepley {
14460598e1a2SLisandro Dalcin   GmshMesh      *mesh = NULL;
144711c56cb0SLisandro Dalcin   PetscViewer    parentviewer = NULL;
14480598e1a2SLisandro Dalcin   PetscBT        periodicVerts = NULL;
14490598e1a2SLisandro Dalcin   PetscBT        periodicCells = NULL;
1450066ea43fSLisandro Dalcin   DM             cdm;
1451331830f3SMatthew G. Knepley   PetscSection   coordSection;
1452331830f3SMatthew G. Knepley   Vec            coordinates;
1453a45dabc8SMatthew G. Knepley   DMLabel        cellSets = NULL, faceSets = NULL, vertSets = NULL, marker = NULL, *regionSets;
1454066ea43fSLisandro Dalcin   PetscInt       dim = 0, coordDim = -1, order = 0;
14550598e1a2SLisandro Dalcin   PetscInt       numNodes = 0, numElems = 0, numVerts = 0, numCells = 0;
1456066ea43fSLisandro Dalcin   PetscInt       cell, cone[8], e, n, v, d;
145781a1af93SMatthew G. Knepley   PetscBool      binary, usemarker = PETSC_FALSE, useregions = PETSC_FALSE, markvertices = PETSC_FALSE;
14580598e1a2SLisandro Dalcin   PetscBool      hybrid = interpolate, periodic = PETSC_TRUE;
1459066ea43fSLisandro Dalcin   PetscBool      highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE;
1460066ea43fSLisandro Dalcin   PetscBool      isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE;
14610598e1a2SLisandro Dalcin   PetscMPIInt    rank;
1462331830f3SMatthew G. Knepley 
1463331830f3SMatthew G. Knepley   PetscFunctionBegin;
14649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
1465d0609cedSBarry Smith   PetscObjectOptionsBegin((PetscObject)viewer);
1466d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"DMPlex Gmsh options");
14679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", hybrid, &hybrid, NULL));
14689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL));
14699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_highorder","Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet));
14709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL));
14719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL));
14729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_regions", "Generate labels with region names", "DMPlexCreateGmsh", useregions, &useregions, NULL));
14739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-dm_plex_gmsh_mark_vertices", "Add vertices to generated labels", "DMPlexCreateGmsh", markvertices, &markvertices, NULL));
14749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE));
1475d0609cedSBarry Smith   PetscOptionsHeadEnd();
1476d0609cedSBarry Smith   PetscOptionsEnd();
14770598e1a2SLisandro Dalcin 
14789566063dSJacob Faibussowitsch   PetscCall(GmshCellInfoSetUp());
147911c56cb0SLisandro Dalcin 
14809566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
14819566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
14829566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL));
148311c56cb0SLisandro Dalcin 
14849566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary));
148511c56cb0SLisandro Dalcin 
148611c56cb0SLisandro Dalcin   /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
14873b46f529SLisandro Dalcin   if (binary) {
148811c56cb0SLisandro Dalcin     parentviewer = viewer;
14899566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer));
149011c56cb0SLisandro Dalcin   }
149111c56cb0SLisandro Dalcin 
1492dd400576SPatrick Sanan   if (rank == 0) {
14930598e1a2SLisandro Dalcin     GmshFile  gmsh[1];
1494698a79b9SLisandro Dalcin     char      line[PETSC_MAX_PATH_LEN];
1495698a79b9SLisandro Dalcin     PetscBool match;
1496331830f3SMatthew G. Knepley 
14979566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(gmsh,1));
1498698a79b9SLisandro Dalcin     gmsh->viewer = viewer;
1499698a79b9SLisandro Dalcin     gmsh->binary = binary;
1500698a79b9SLisandro Dalcin 
15019566063dSJacob Faibussowitsch     PetscCall(GmshMeshCreate(&mesh));
15020598e1a2SLisandro Dalcin 
1503698a79b9SLisandro Dalcin     /* Read mesh format */
15049566063dSJacob Faibussowitsch     PetscCall(GmshReadSection(gmsh, line));
15059566063dSJacob Faibussowitsch     PetscCall(GmshExpect(gmsh, "$MeshFormat", line));
15069566063dSJacob Faibussowitsch     PetscCall(GmshReadMeshFormat(gmsh));
15079566063dSJacob Faibussowitsch     PetscCall(GmshReadEndSection(gmsh, "$EndMeshFormat", line));
150811c56cb0SLisandro Dalcin 
1509dc0ede3bSMatthew G. Knepley     /* OPTIONAL Read physical names */
15109566063dSJacob Faibussowitsch     PetscCall(GmshReadSection(gmsh, line));
15119566063dSJacob Faibussowitsch     PetscCall(GmshMatch(gmsh, "$PhysicalNames", line, &match));
1512dc0ede3bSMatthew G. Knepley     if (match) {
15139566063dSJacob Faibussowitsch       PetscCall(GmshExpect(gmsh, "$PhysicalNames", line));
15149566063dSJacob Faibussowitsch       PetscCall(GmshReadPhysicalNames(gmsh, mesh));
15159566063dSJacob Faibussowitsch       PetscCall(GmshReadEndSection(gmsh, "$EndPhysicalNames", line));
1516698a79b9SLisandro Dalcin       /* Initial read for entity section */
15179566063dSJacob Faibussowitsch       PetscCall(GmshReadSection(gmsh, line));
1518dc0ede3bSMatthew G. Knepley     }
151911c56cb0SLisandro Dalcin 
1520de78e4feSLisandro Dalcin     /* Read entities */
1521698a79b9SLisandro Dalcin     if (gmsh->fileFormat >= 40) {
15229566063dSJacob Faibussowitsch       PetscCall(GmshExpect(gmsh, "$Entities", line));
15239566063dSJacob Faibussowitsch       PetscCall(GmshReadEntities(gmsh, mesh));
15249566063dSJacob Faibussowitsch       PetscCall(GmshReadEndSection(gmsh, "$EndEntities", line));
1525698a79b9SLisandro Dalcin       /* Initial read for nodes section */
15269566063dSJacob Faibussowitsch       PetscCall(GmshReadSection(gmsh, line));
1527de78e4feSLisandro Dalcin     }
1528de78e4feSLisandro Dalcin 
1529698a79b9SLisandro Dalcin     /* Read nodes */
15309566063dSJacob Faibussowitsch     PetscCall(GmshExpect(gmsh, "$Nodes", line));
15319566063dSJacob Faibussowitsch     PetscCall(GmshReadNodes(gmsh, mesh));
15329566063dSJacob Faibussowitsch     PetscCall(GmshReadEndSection(gmsh, "$EndNodes", line));
153311c56cb0SLisandro Dalcin 
1534698a79b9SLisandro Dalcin     /* Read elements */
15359566063dSJacob Faibussowitsch     PetscCall(GmshReadSection(gmsh, line));
15369566063dSJacob Faibussowitsch     PetscCall(GmshExpect(gmsh, "$Elements", line));
15379566063dSJacob Faibussowitsch     PetscCall(GmshReadElements(gmsh, mesh));
15389566063dSJacob Faibussowitsch     PetscCall(GmshReadEndSection(gmsh, "$EndElements", line));
1539de78e4feSLisandro Dalcin 
15400598e1a2SLisandro Dalcin     /* Read periodic section (OPTIONAL) */
1541698a79b9SLisandro Dalcin     if (periodic) {
15429566063dSJacob Faibussowitsch       PetscCall(GmshReadSection(gmsh, line));
15439566063dSJacob Faibussowitsch       PetscCall(GmshMatch(gmsh, "$Periodic", line, &periodic));
1544ef12879bSLisandro Dalcin     }
1545ef12879bSLisandro Dalcin     if (periodic) {
15469566063dSJacob Faibussowitsch       PetscCall(GmshExpect(gmsh, "$Periodic", line));
15479566063dSJacob Faibussowitsch       PetscCall(GmshReadPeriodic(gmsh, mesh));
15489566063dSJacob Faibussowitsch       PetscCall(GmshReadEndSection(gmsh, "$EndPeriodic", line));
1549698a79b9SLisandro Dalcin     }
1550698a79b9SLisandro Dalcin 
15519566063dSJacob Faibussowitsch     PetscCall(PetscFree(gmsh->wbuf));
15529566063dSJacob Faibussowitsch     PetscCall(PetscFree(gmsh->sbuf));
15539566063dSJacob Faibussowitsch     PetscCall(PetscFree(gmsh->nbuf));
15540598e1a2SLisandro Dalcin 
15550598e1a2SLisandro Dalcin     dim       = mesh->dim;
1556066ea43fSLisandro Dalcin     order     = mesh->order;
15570598e1a2SLisandro Dalcin     numNodes  = mesh->numNodes;
15580598e1a2SLisandro Dalcin     numElems  = mesh->numElems;
15590598e1a2SLisandro Dalcin     numVerts  = mesh->numVerts;
15600598e1a2SLisandro Dalcin     numCells  = mesh->numCells;
1561066ea43fSLisandro Dalcin 
1562066ea43fSLisandro Dalcin     {
1563066ea43fSLisandro Dalcin       GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL;
1564066ea43fSLisandro Dalcin       GmshElement *elemB = elemA ? elemA + mesh->numCells - 1  : NULL;
1565066ea43fSLisandro Dalcin       int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1;
1566066ea43fSLisandro Dalcin       int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1;
1567066ea43fSLisandro Dalcin       isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE;
1568066ea43fSLisandro Dalcin       isHybrid  = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE;
1569066ea43fSLisandro Dalcin       hasTetra  = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE;
1570066ea43fSLisandro Dalcin     }
1571698a79b9SLisandro Dalcin   }
1572698a79b9SLisandro Dalcin 
1573698a79b9SLisandro Dalcin   if (parentviewer) {
15749566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer));
1575698a79b9SLisandro Dalcin   }
1576698a79b9SLisandro Dalcin 
1577066ea43fSLisandro Dalcin   {
1578066ea43fSLisandro Dalcin     int buf[6];
1579698a79b9SLisandro Dalcin 
1580066ea43fSLisandro Dalcin     buf[0] = (int)dim;
1581066ea43fSLisandro Dalcin     buf[1] = (int)order;
1582066ea43fSLisandro Dalcin     buf[2] = periodic;
1583066ea43fSLisandro Dalcin     buf[3] = isSimplex;
1584066ea43fSLisandro Dalcin     buf[4] = isHybrid;
1585066ea43fSLisandro Dalcin     buf[5] = hasTetra;
1586066ea43fSLisandro Dalcin 
15879566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(buf, 6, MPI_INT, 0, comm));
1588066ea43fSLisandro Dalcin 
1589066ea43fSLisandro Dalcin     dim       = buf[0];
1590066ea43fSLisandro Dalcin     order     = buf[1];
1591066ea43fSLisandro Dalcin     periodic  = buf[2] ? PETSC_TRUE : PETSC_FALSE;
1592066ea43fSLisandro Dalcin     isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE;
1593066ea43fSLisandro Dalcin     isHybrid  = buf[4] ? PETSC_TRUE : PETSC_FALSE;
1594066ea43fSLisandro Dalcin     hasTetra  = buf[5] ? PETSC_TRUE : PETSC_FALSE;
1595dea2a3c7SStefano Zampini   }
1596dea2a3c7SStefano Zampini 
1597066ea43fSLisandro Dalcin   if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE;
159808401ef6SPierre Jolivet   PetscCheck(!highOrder || !isHybrid,comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet");
1599066ea43fSLisandro Dalcin 
16000598e1a2SLisandro Dalcin   /* We do not want this label automatically computed, instead we fill it here */
16019566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "celltype"));
160211c56cb0SLisandro Dalcin 
1603a4bb7517SMichael Lange   /* Allocate the cell-vertex mesh */
16049566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, numCells+numVerts));
16050598e1a2SLisandro Dalcin   for (cell = 0; cell < numCells; ++cell) {
16060598e1a2SLisandro Dalcin     GmshElement *elem = mesh->elements + cell;
16070598e1a2SLisandro Dalcin     DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType);
16080598e1a2SLisandro Dalcin     if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR;
16099566063dSJacob Faibussowitsch     PetscCall(DMPlexSetConeSize(*dm, cell, elem->numVerts));
16109566063dSJacob Faibussowitsch     PetscCall(DMPlexSetCellType(*dm, cell, ctype));
1611db397164SMichael Lange   }
16120598e1a2SLisandro Dalcin   for (v = numCells; v < numCells+numVerts; ++v) {
16139566063dSJacob Faibussowitsch     PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
161496ca5757SLisandro Dalcin   }
16159566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*dm));
161696ca5757SLisandro Dalcin 
1617a4bb7517SMichael Lange   /* Add cell-vertex connections */
16180598e1a2SLisandro Dalcin   for (cell = 0; cell < numCells; ++cell) {
16190598e1a2SLisandro Dalcin     GmshElement *elem = mesh->elements + cell;
16200598e1a2SLisandro Dalcin     for (v = 0; v < elem->numVerts; ++v) {
16210598e1a2SLisandro Dalcin       const PetscInt nn = elem->nodes[v];
16220598e1a2SLisandro Dalcin       const PetscInt vv = mesh->vertexMap[nn];
16230598e1a2SLisandro Dalcin       cone[v] = numCells + vv;
1624db397164SMichael Lange     }
16259566063dSJacob Faibussowitsch     PetscCall(DMPlexReorderCell(*dm, cell, cone));
16269566063dSJacob Faibussowitsch     PetscCall(DMPlexSetCone(*dm, cell, cone));
1627a4bb7517SMichael Lange   }
162896ca5757SLisandro Dalcin 
16299566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
16309566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
16319566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
1632331830f3SMatthew G. Knepley   if (interpolate) {
16335fd9971aSMatthew G. Knepley     DM idm;
1634331830f3SMatthew G. Knepley 
16359566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
16369566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1637331830f3SMatthew G. Knepley     *dm  = idm;
1638331830f3SMatthew G. Knepley   }
1639917266a4SLisandro Dalcin 
164081a1af93SMatthew G. Knepley   /* Create the label "marker" over the whole boundary */
16412c71b3e2SJacob Faibussowitsch   PetscCheckFalse(usemarker && !interpolate && dim > 1,comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation");
1642dd400576SPatrick Sanan   if (rank == 0 && usemarker) {
1643d3f73514SStefano Zampini     PetscInt f, fStart, fEnd;
1644d3f73514SStefano Zampini 
16459566063dSJacob Faibussowitsch     PetscCall(DMCreateLabel(*dm, "marker"));
16469566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd));
1647d3f73514SStefano Zampini     for (f = fStart; f < fEnd; ++f) {
1648d3f73514SStefano Zampini       PetscInt suppSize;
1649d3f73514SStefano Zampini 
16509566063dSJacob Faibussowitsch       PetscCall(DMPlexGetSupportSize(*dm, f, &suppSize));
1651d3f73514SStefano Zampini       if (suppSize == 1) {
1652d3f73514SStefano Zampini         PetscInt *cone = NULL, coneSize, p;
1653d3f73514SStefano Zampini 
16549566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone));
1655d3f73514SStefano Zampini         for (p = 0; p < coneSize; p += 2) {
16569566063dSJacob Faibussowitsch           PetscCall(DMSetLabelValue_Fast(*dm, &marker, "marker", cone[p], 1));
1657d3f73514SStefano Zampini         }
16589566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone));
1659d3f73514SStefano Zampini       }
1660d3f73514SStefano Zampini     }
1661d3f73514SStefano Zampini   }
166216ad7e67SMichael Lange 
1663dd400576SPatrick Sanan   if (rank == 0) {
1664a45dabc8SMatthew G. Knepley     const PetscInt Nr = useregions ? mesh->numRegions : 0;
166511c56cb0SLisandro Dalcin     PetscInt       vStart, vEnd;
1666d1a54cefSMatthew G. Knepley 
16679566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(Nr, &regionSets));
16689566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd));
16690598e1a2SLisandro Dalcin     for (cell = 0, e = 0; e < numElems; ++e) {
16700598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + e;
16716e1dd89aSLawrence Mitchell 
16726e1dd89aSLawrence Mitchell       /* Create cell sets */
16730598e1a2SLisandro Dalcin       if (elem->dim == dim && dim > 0) {
16740598e1a2SLisandro Dalcin         if (elem->numTags > 0) {
1675a45dabc8SMatthew G. Knepley           const PetscInt tag = elem->tags[0];
1676a45dabc8SMatthew G. Knepley           PetscInt       r;
1677a45dabc8SMatthew G. Knepley 
16789566063dSJacob Faibussowitsch           if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", cell, tag));
1679a45dabc8SMatthew G. Knepley           for (r = 0; r < Nr; ++r) {
16809566063dSJacob Faibussowitsch             if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, &regionSets[r], mesh->regionNames[r], cell, tag));
1681a45dabc8SMatthew G. Knepley           }
16826e1dd89aSLawrence Mitchell         }
1683917266a4SLisandro Dalcin         cell++;
16846e1dd89aSLawrence Mitchell       }
16850c070f12SSander Arens 
16860598e1a2SLisandro Dalcin       /* Create face sets */
16870598e1a2SLisandro Dalcin       if (interpolate && elem->dim == dim-1) {
16880598e1a2SLisandro Dalcin         PetscInt        joinSize;
16890598e1a2SLisandro Dalcin         const PetscInt *join = NULL;
1690a45dabc8SMatthew G. Knepley         const PetscInt  tag = elem->tags[0];
1691a45dabc8SMatthew G. Knepley         PetscInt        r;
1692a45dabc8SMatthew G. Knepley 
16930598e1a2SLisandro Dalcin         /* Find the relevant facet with vertex joins */
16940598e1a2SLisandro Dalcin         for (v = 0; v < elem->numVerts; ++v) {
16950598e1a2SLisandro Dalcin           const PetscInt nn = elem->nodes[v];
16960598e1a2SLisandro Dalcin           const PetscInt vv = mesh->vertexMap[nn];
16970598e1a2SLisandro Dalcin           cone[v] = vStart + vv;
16980598e1a2SLisandro Dalcin         }
16999566063dSJacob Faibussowitsch         PetscCall(DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join));
170008401ef6SPierre Jolivet         PetscCheck(joinSize == 1,PETSC_COMM_SELF, PETSC_ERR_SUP, "Could not determine Plex facet for Gmsh element %D (Plex cell %D)", elem->id, e);
17019566063dSJacob Faibussowitsch         if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], tag));
1702a45dabc8SMatthew G. Knepley         for (r = 0; r < Nr; ++r) {
17039566063dSJacob Faibussowitsch           if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, &regionSets[r], mesh->regionNames[r], join[0], tag));
1704a45dabc8SMatthew G. Knepley         }
17059566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join));
17060598e1a2SLisandro Dalcin       }
17070598e1a2SLisandro Dalcin 
17080c070f12SSander Arens       /* Create vertex sets */
17090598e1a2SLisandro Dalcin       if (elem->dim == 0) {
17100598e1a2SLisandro Dalcin         if (elem->numTags > 0) {
17110598e1a2SLisandro Dalcin           const PetscInt nn = elem->nodes[0];
17120598e1a2SLisandro Dalcin           const PetscInt vv = mesh->vertexMap[nn];
1713a45dabc8SMatthew G. Knepley           const PetscInt tag = elem->tags[0];
1714a45dabc8SMatthew G. Knepley           PetscInt       r;
1715a45dabc8SMatthew G. Knepley 
17169566063dSJacob Faibussowitsch           if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag));
171781a1af93SMatthew G. Knepley           for (r = 0; r < Nr; ++r) {
17189566063dSJacob Faibussowitsch             if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, &regionSets[r], mesh->regionNames[r], vStart + vv, tag));
171981a1af93SMatthew G. Knepley           }
172081a1af93SMatthew G. Knepley         }
172181a1af93SMatthew G. Knepley       }
172281a1af93SMatthew G. Knepley     }
172381a1af93SMatthew G. Knepley     if (markvertices) {
172481a1af93SMatthew G. Knepley       for (v = 0; v < numNodes; ++v) {
172581a1af93SMatthew G. Knepley         const PetscInt vv  = mesh->vertexMap[v];
172681a1af93SMatthew G. Knepley         const PetscInt tag = mesh->nodelist->tag[v];
172781a1af93SMatthew G. Knepley         PetscInt       r;
172881a1af93SMatthew G. Knepley 
172981a1af93SMatthew G. Knepley         if (tag != -1) {
17309566063dSJacob Faibussowitsch           if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag));
1731a45dabc8SMatthew G. Knepley           for (r = 0; r < Nr; ++r) {
17329566063dSJacob Faibussowitsch             if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, &regionSets[r], mesh->regionNames[r], vStart + vv, tag));
17330598e1a2SLisandro Dalcin           }
17340598e1a2SLisandro Dalcin         }
17350598e1a2SLisandro Dalcin       }
17360598e1a2SLisandro Dalcin     }
17379566063dSJacob Faibussowitsch     PetscCall(PetscFree(regionSets));
1738a45dabc8SMatthew G. Knepley   }
17390598e1a2SLisandro Dalcin 
17407dd454faSLisandro Dalcin   { /* Create Cell/Face/Vertex Sets labels at all processes */
17417dd454faSLisandro Dalcin     enum {n = 4};
17427dd454faSLisandro Dalcin     PetscBool flag[n];
17437dd454faSLisandro Dalcin 
17447dd454faSLisandro Dalcin     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
17457dd454faSLisandro Dalcin     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
17467dd454faSLisandro Dalcin     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
17477dd454faSLisandro Dalcin     flag[3] = marker   ? PETSC_TRUE : PETSC_FALSE;
17489566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
17499566063dSJacob Faibussowitsch     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
17509566063dSJacob Faibussowitsch     if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
17519566063dSJacob Faibussowitsch     if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
17529566063dSJacob Faibussowitsch     if (flag[3]) PetscCall(DMCreateLabel(*dm, "marker"));
17537dd454faSLisandro Dalcin   }
17547dd454faSLisandro Dalcin 
17550598e1a2SLisandro Dalcin   if (periodic) {
17569566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(numVerts, &periodicVerts));
17570598e1a2SLisandro Dalcin     for (n = 0; n < numNodes; ++n) {
17580598e1a2SLisandro Dalcin       if (mesh->vertexMap[n] >= 0) {
17590598e1a2SLisandro Dalcin         if (PetscUnlikely(mesh->periodMap[n] != n)) {
17600598e1a2SLisandro Dalcin           PetscInt m = mesh->periodMap[n];
17619566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[n]));
17629566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[m]));
17630598e1a2SLisandro Dalcin         }
17640598e1a2SLisandro Dalcin       }
17650598e1a2SLisandro Dalcin     }
17669566063dSJacob Faibussowitsch     PetscCall(PetscBTCreate(numCells, &periodicCells));
17670598e1a2SLisandro Dalcin     for (cell = 0; cell < numCells; ++cell) {
17680598e1a2SLisandro Dalcin       GmshElement *elem = mesh->elements + cell;
17690598e1a2SLisandro Dalcin       for (v = 0; v < elem->numVerts; ++v) {
17700598e1a2SLisandro Dalcin         PetscInt nn = elem->nodes[v];
17710598e1a2SLisandro Dalcin         PetscInt vv = mesh->vertexMap[nn];
17720598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) {
17739566063dSJacob Faibussowitsch           PetscCall(PetscBTSet(periodicCells, cell)); break;
17740c070f12SSander Arens         }
17750c070f12SSander Arens       }
17760c070f12SSander Arens     }
177716ad7e67SMichael Lange   }
177816ad7e67SMichael Lange 
1779066ea43fSLisandro Dalcin   /* Setup coordinate DM */
17800598e1a2SLisandro Dalcin   if (coordDim < 0) coordDim = dim;
17819566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDim(*dm, coordDim));
17829566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*dm, &cdm));
1783066ea43fSLisandro Dalcin   if (highOrder) {
1784066ea43fSLisandro Dalcin     PetscFE         fe;
1785066ea43fSLisandro Dalcin     PetscBool       continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
1786066ea43fSLisandro Dalcin     PetscDTNodeType nodeType   = PETSCDTNODES_EQUISPACED;
1787066ea43fSLisandro Dalcin 
1788066ea43fSLisandro Dalcin     if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */
1789066ea43fSLisandro Dalcin 
17909566063dSJacob Faibussowitsch     PetscCall(GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe));
17919566063dSJacob Faibussowitsch     PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view"));
17929566063dSJacob Faibussowitsch     PetscCall(DMSetField(cdm, 0, NULL, (PetscObject) fe));
17939566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
17949566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(cdm));
1795066ea43fSLisandro Dalcin   }
1796066ea43fSLisandro Dalcin 
1797066ea43fSLisandro Dalcin   /* Create coordinates */
1798066ea43fSLisandro Dalcin   if (highOrder) {
1799066ea43fSLisandro Dalcin 
1800066ea43fSLisandro Dalcin     PetscInt     maxDof = GmshNumNodes_HEX(order)*coordDim;
1801066ea43fSLisandro Dalcin     double       *coords = mesh ? mesh->nodelist->xyz : NULL;
1802066ea43fSLisandro Dalcin     PetscSection section;
1803066ea43fSLisandro Dalcin     PetscScalar  *cellCoords;
1804066ea43fSLisandro Dalcin 
18059566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(cdm, NULL));
18069566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(cdm, &coordSection));
18079566063dSJacob Faibussowitsch     PetscCall(PetscSectionClone(coordSection, &section));
18089566063dSJacob Faibussowitsch     PetscCall(DMPlexSetClosurePermutationTensor(cdm, 0, section)); /* XXX Implement DMPlexSetClosurePermutationLexicographic() */
1809066ea43fSLisandro Dalcin 
18109566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(cdm, &coordinates));
18119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(maxDof, &cellCoords));
1812066ea43fSLisandro Dalcin     for (cell = 0; cell < numCells; ++cell) {
1813066ea43fSLisandro Dalcin       GmshElement *elem = mesh->elements + cell;
1814066ea43fSLisandro Dalcin       const int *lexorder = GmshCellMap[elem->cellType].lexorder();
1815b9bf55e5SMatthew G. Knepley       int s = 0;
1816066ea43fSLisandro Dalcin       for (n = 0; n < elem->numNodes; ++n) {
1817b9bf55e5SMatthew G. Knepley         while (lexorder[n+s] < 0) ++s;
1818b9bf55e5SMatthew G. Knepley         const PetscInt node = elem->nodes[lexorder[n+s]];
1819b9bf55e5SMatthew G. Knepley         for (d = 0; d < coordDim; ++d) cellCoords[(n+s)*coordDim+d] = (PetscReal) coords[node*3+d];
1820b9bf55e5SMatthew G. Knepley       }
1821b9bf55e5SMatthew G. Knepley       if (s) {
1822b9bf55e5SMatthew G. Knepley         /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */
1823b9bf55e5SMatthew G. Knepley         PetscReal quaCenterWeights[9]  = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25};
1824b9bf55e5SMatthew G. Knepley         /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */
1825b9bf55e5SMatthew G. Knepley         PetscReal hexBottomWeights[27] = {-0.25, 0.5,  -0.25, 0.5,  0.0, 0.5,  -0.25, 0.5,  -0.25,
1826b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.0,  0.0,  0.0, 0.0,   0.0,  0.0,   0.0,
1827b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.0,  0.0,  0.0, 0.0,   0.0,  0.0,   0.0};
1828b9bf55e5SMatthew G. Knepley         PetscReal hexFrontWeights[27]  = {-0.25, 0.5,  -0.25, 0.0,  0.0, 0.0,   0.0,  0.0,   0.0,
1829b9bf55e5SMatthew G. Knepley                                            0.5,  0.0,   0.5,  0.0,  0.0, 0.0,   0.0,  0.0,   0.0,
1830b9bf55e5SMatthew G. Knepley                                           -0.25, 0.5,  -0.25, 0.0,  0.0, 0.0,   0.0,  0.0,   0.0};
1831b9bf55e5SMatthew G. Knepley         PetscReal hexLeftWeights[27]   = {-0.25, 0.0,   0.0,  0.5,  0.0, 0.0,  -0.25, 0.0,   0.0,
1832b9bf55e5SMatthew G. Knepley                                            0.5,  0.0,   0.0,  0.0,  0.0, 0.0,   0.5,  0.0,   0.0,
1833b9bf55e5SMatthew G. Knepley                                           -0.25, 0.0,   0.0,  0.5,  0.0, 0.0,  -0.25, 0.0,   0.0};
1834b9bf55e5SMatthew G. Knepley         PetscReal hexRightWeights[27]  = { 0.0,  0.0,  -0.25, 0.0,  0.0, 0.5,   0.0,  0.0,  -0.25,
1835b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.5,  0.0,  0.0, 0.0,   0.0,  0.0,   0.5,
1836b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,  -0.25, 0.0,  0.0, 0.5,   0.0,  0.0,  -0.25};
1837b9bf55e5SMatthew G. Knepley         PetscReal hexBackWeights[27]   = { 0.0,  0.0,   0.0,  0.0,  0.0, 0.0,  -0.25, 0.5,  -0.25,
1838b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.0,  0.0,  0.0, 0.0,   0.5,  0.0,   0.5,
1839b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.0,  0.0,  0.0, 0.0,  -0.25, 0.5,  -0.25};
1840b9bf55e5SMatthew G. Knepley         PetscReal hexTopWeights[27]    = { 0.0,  0.0,   0.0,  0.0,  0.0, 0.0,   0.0,  0.0,   0.0,
1841b9bf55e5SMatthew G. Knepley                                            0.0,  0.0,   0.0,  0.0,  0.0, 0.0,   0.0,  0.0,   0.0,
1842b9bf55e5SMatthew G. Knepley                                           -0.25, 0.5,  -0.25, 0.5,  0.0, 0.5,  -0.25, 0.5,  -0.25};
1843b9bf55e5SMatthew G. Knepley         PetscReal hexCenterWeights[27] = {-0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25,
1844b9bf55e5SMatthew G. Knepley                                            0.25, 0.0,   0.25, 0.0,  0.0, 0.0,   0.25, 0.0,   0.25,
1845b9bf55e5SMatthew G. Knepley                                           -0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25};
1846b9bf55e5SMatthew G. Knepley         PetscReal  *sdWeights2[9]      = {NULL, NULL, NULL, NULL, quaCenterWeights, NULL, NULL, NULL, NULL};
1847b9bf55e5SMatthew G. Knepley         PetscReal  *sdWeights3[27]     = {NULL, NULL, NULL, NULL, hexBottomWeights, NULL, NULL, NULL, NULL,
1848b9bf55e5SMatthew G. Knepley                                           NULL, hexFrontWeights, NULL, hexLeftWeights, hexCenterWeights, hexRightWeights, NULL, hexBackWeights, NULL,
1849b9bf55e5SMatthew G. Knepley                                           NULL, NULL, NULL, NULL, hexTopWeights,    NULL, NULL, NULL, NULL};
1850b9bf55e5SMatthew G. Knepley         PetscReal **sdWeights[4]       = {NULL, NULL, sdWeights2, sdWeights3};
1851b9bf55e5SMatthew G. Knepley 
1852b9bf55e5SMatthew G. Knepley         /* Missing entries in serendipity cell, only works for 8-node quad and 20-node hex */
1853b9bf55e5SMatthew G. Knepley         for (n = 0; n < elem->numNodes+s; ++n) {
1854b9bf55e5SMatthew G. Knepley           if (lexorder[n] >= 0) continue;
1855b9bf55e5SMatthew G. Knepley           for (d = 0; d < coordDim; ++d) cellCoords[n*coordDim+d] = 0.0;
1856b9bf55e5SMatthew G. Knepley           for (int bn = 0; bn < elem->numNodes+s; ++bn) {
1857b9bf55e5SMatthew G. Knepley             if (lexorder[bn] < 0) continue;
1858b9bf55e5SMatthew G. Knepley             const PetscReal *weights = sdWeights[coordDim][n];
1859b9bf55e5SMatthew G. Knepley             const PetscInt   bnode   = elem->nodes[lexorder[bn]];
1860b9bf55e5SMatthew G. Knepley             for (d = 0; d < coordDim; ++d) cellCoords[n*coordDim+d] += weights[bn] * (PetscReal) coords[bnode*3+d];
1861b9bf55e5SMatthew G. Knepley           }
1862b9bf55e5SMatthew G. Knepley         }
1863066ea43fSLisandro Dalcin       }
18649566063dSJacob Faibussowitsch       PetscCall(DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES));
1865066ea43fSLisandro Dalcin     }
18669566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&section));
18679566063dSJacob Faibussowitsch     PetscCall(PetscFree(cellCoords));
1868066ea43fSLisandro Dalcin 
1869066ea43fSLisandro Dalcin   } else {
1870066ea43fSLisandro Dalcin 
1871066ea43fSLisandro Dalcin     PetscInt    *nodeMap;
1872066ea43fSLisandro Dalcin     double      *coords = mesh ? mesh->nodelist->xyz : NULL;
1873066ea43fSLisandro Dalcin     PetscScalar *pointCoords;
1874066ea43fSLisandro Dalcin 
18759566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(cdm, &coordSection));
18769566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetNumFields(coordSection, 1));
18779566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
18780598e1a2SLisandro Dalcin     if (periodic) { /* we need to localize coordinates on cells */
18799566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetChart(coordSection, 0, numCells+numVerts));
1880f45c9589SStefano Zampini     } else {
18819566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetChart(coordSection, numCells, numCells+numVerts));
1882f45c9589SStefano Zampini     }
18830598e1a2SLisandro Dalcin     if (periodic) {
18840598e1a2SLisandro Dalcin       for (cell = 0; cell < numCells; ++cell) {
18850598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) {
18860598e1a2SLisandro Dalcin           GmshElement *elem = mesh->elements + cell;
18870598e1a2SLisandro Dalcin           PetscInt dof = elem->numVerts * coordDim;
18889566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetDof(coordSection, cell, dof));
18899566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetFieldDof(coordSection, cell, 0, dof));
1890f45c9589SStefano Zampini         }
1891f45c9589SStefano Zampini       }
1892f45c9589SStefano Zampini     }
18930598e1a2SLisandro Dalcin     for (v = numCells; v < numCells+numVerts; ++v) {
18949566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(coordSection, v, coordDim));
18959566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
18960598e1a2SLisandro Dalcin     }
18979566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(coordSection));
1898066ea43fSLisandro Dalcin 
18999566063dSJacob Faibussowitsch     PetscCall(DMCreateLocalVector(cdm, &coordinates));
19009566063dSJacob Faibussowitsch     PetscCall(VecGetArray(coordinates, &pointCoords));
19010598e1a2SLisandro Dalcin     if (periodic) {
19020598e1a2SLisandro Dalcin       for (cell = 0; cell < numCells; ++cell) {
19030598e1a2SLisandro Dalcin         if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) {
19040598e1a2SLisandro Dalcin           GmshElement *elem = mesh->elements + cell;
19050598e1a2SLisandro Dalcin           PetscInt off, node;
19060598e1a2SLisandro Dalcin           for (v = 0; v < elem->numVerts; ++v)
19070598e1a2SLisandro Dalcin             cone[v] = elem->nodes[v];
19089566063dSJacob Faibussowitsch           PetscCall(DMPlexReorderCell(cdm, cell, cone));
19099566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetOffset(coordSection, cell, &off));
19100598e1a2SLisandro Dalcin           for (v = 0; v < elem->numVerts; ++v)
19110598e1a2SLisandro Dalcin             for (node = cone[v], d = 0; d < coordDim; ++d)
1912066ea43fSLisandro Dalcin               pointCoords[off++] = (PetscReal) coords[node*3+d];
1913331830f3SMatthew G. Knepley         }
1914331830f3SMatthew G. Knepley       }
1915331830f3SMatthew G. Knepley     }
19169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numVerts, &nodeMap));
19170598e1a2SLisandro Dalcin     for (n = 0; n < numNodes; n++)
19180598e1a2SLisandro Dalcin       if (mesh->vertexMap[n] >= 0)
1919066ea43fSLisandro Dalcin         nodeMap[mesh->vertexMap[n]] = n;
19200598e1a2SLisandro Dalcin     for (v = 0; v < numVerts; ++v) {
1921066ea43fSLisandro Dalcin       PetscInt off, node = nodeMap[v];
19229566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(coordSection, numCells + v, &off));
19230598e1a2SLisandro Dalcin       for (d = 0; d < coordDim; ++d)
1924066ea43fSLisandro Dalcin         pointCoords[off+d] = (PetscReal) coords[node*3+d];
19250598e1a2SLisandro Dalcin     }
19269566063dSJacob Faibussowitsch     PetscCall(PetscFree(nodeMap));
19279566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(coordinates, &pointCoords));
1928066ea43fSLisandro Dalcin 
1929066ea43fSLisandro Dalcin   }
1930066ea43fSLisandro Dalcin 
19319566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
19329566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, coordDim));
19339566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
19349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
19359566063dSJacob Faibussowitsch   PetscCall(DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL));
193611c56cb0SLisandro Dalcin 
19379566063dSJacob Faibussowitsch   PetscCall(GmshMeshDestroy(&mesh));
19389566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&periodicVerts));
19399566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&periodicCells));
194011c56cb0SLisandro Dalcin 
1941066ea43fSLisandro Dalcin   if (highOrder && project)  {
1942066ea43fSLisandro Dalcin     PetscFE         fe;
1943066ea43fSLisandro Dalcin     const char      prefix[]   = "dm_plex_gmsh_project_";
1944066ea43fSLisandro Dalcin     PetscBool       continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
1945066ea43fSLisandro Dalcin     PetscDTNodeType nodeType   = PETSCDTNODES_GAUSSJACOBI;
1946066ea43fSLisandro Dalcin 
1947066ea43fSLisandro Dalcin     if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */
1948066ea43fSLisandro Dalcin 
19499566063dSJacob Faibussowitsch     PetscCall(GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe));
19509566063dSJacob Faibussowitsch     PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view"));
19519566063dSJacob Faibussowitsch     PetscCall(DMProjectCoordinates(*dm, fe));
19529566063dSJacob Faibussowitsch     PetscCall(PetscFEDestroy(&fe));
1953066ea43fSLisandro Dalcin   }
1954066ea43fSLisandro Dalcin 
19559566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL));
1956331830f3SMatthew G. Knepley   PetscFunctionReturn(0);
1957331830f3SMatthew G. Knepley }
1958