Actual source code: plexgmsh.c
1: #include <petsc/private/dmpleximpl.h>
2: #include <petsc/private/hashmapi.h>
4: #include <../src/dm/impls/plex/gmshlex.h>
6: #define GMSH_LEXORDER_ITEM(T, p) \
7: static int *GmshLexOrder_##T##_##p(void) \
8: { \
9: static int Gmsh_LexOrder_##T##_##p[GmshNumNodes_##T(p)] = {-1}; \
10: int *lex = Gmsh_LexOrder_##T##_##p; \
11: if (lex[0] == -1) (void)GmshLexOrder_##T(p, lex, 0); \
12: return lex; \
13: }
15: static int *GmshLexOrder_QUA_2_Serendipity(void)
16: {
17: static int Gmsh_LexOrder_QUA_2_Serendipity[9] = {-1};
18: int *lex = Gmsh_LexOrder_QUA_2_Serendipity;
19: if (lex[0] == -1) {
20: /* Vertices */
21: lex[0] = 0;
22: lex[2] = 1;
23: lex[8] = 2;
24: lex[6] = 3;
25: /* Edges */
26: lex[1] = 4;
27: lex[5] = 5;
28: lex[7] = 6;
29: lex[3] = 7;
30: /* Cell */
31: lex[4] = -8;
32: }
33: return lex;
34: }
36: static int *GmshLexOrder_HEX_2_Serendipity(void)
37: {
38: static int Gmsh_LexOrder_HEX_2_Serendipity[27] = {-1};
39: int *lex = Gmsh_LexOrder_HEX_2_Serendipity;
40: if (lex[0] == -1) {
41: /* Vertices */
42: lex[0] = 0;
43: lex[2] = 1;
44: lex[8] = 2;
45: lex[6] = 3;
46: lex[18] = 4;
47: lex[20] = 5;
48: lex[26] = 6;
49: lex[24] = 7;
50: /* Edges */
51: lex[1] = 8;
52: lex[3] = 9;
53: lex[9] = 10;
54: lex[5] = 11;
55: lex[11] = 12;
56: lex[7] = 13;
57: lex[17] = 14;
58: lex[15] = 15;
59: lex[19] = 16;
60: lex[21] = 17;
61: lex[23] = 18;
62: lex[25] = 19;
63: /* Faces */
64: lex[4] = -20;
65: lex[10] = -21;
66: lex[12] = -22;
67: lex[14] = -23;
68: lex[16] = -24;
69: lex[22] = -25;
70: /* Cell */
71: lex[13] = -26;
72: }
73: return lex;
74: }
76: #define GMSH_LEXORDER_LIST(T) \
77: GMSH_LEXORDER_ITEM(T, 1) \
78: GMSH_LEXORDER_ITEM(T, 2) \
79: GMSH_LEXORDER_ITEM(T, 3) \
80: GMSH_LEXORDER_ITEM(T, 4) \
81: GMSH_LEXORDER_ITEM(T, 5) \
82: GMSH_LEXORDER_ITEM(T, 6) \
83: GMSH_LEXORDER_ITEM(T, 7) \
84: GMSH_LEXORDER_ITEM(T, 8) \
85: GMSH_LEXORDER_ITEM(T, 9) \
86: GMSH_LEXORDER_ITEM(T, 10)
88: GMSH_LEXORDER_ITEM(VTX, 0)
89: GMSH_LEXORDER_LIST(SEG)
90: GMSH_LEXORDER_LIST(TRI)
91: GMSH_LEXORDER_LIST(QUA)
92: GMSH_LEXORDER_LIST(TET)
93: GMSH_LEXORDER_LIST(HEX)
94: GMSH_LEXORDER_LIST(PRI)
95: GMSH_LEXORDER_LIST(PYR)
97: typedef enum {
98: GMSH_VTX = 0,
99: GMSH_SEG = 1,
100: GMSH_TRI = 2,
101: GMSH_QUA = 3,
102: GMSH_TET = 4,
103: GMSH_HEX = 5,
104: GMSH_PRI = 6,
105: GMSH_PYR = 7,
106: GMSH_NUM_POLYTOPES = 8
107: } GmshPolytopeType;
109: typedef struct {
110: int cellType;
111: int polytope;
112: int dim;
113: int order;
114: int numVerts;
115: int numNodes;
116: int *(*lexorder)(void);
117: } GmshCellInfo;
119: #define GmshCellEntry(cellType, polytope, dim, order) {cellType, GMSH_##polytope, dim, order, GmshNumNodes_##polytope(1), GmshNumNodes_##polytope(order), GmshLexOrder_##polytope##_##order}
121: static const GmshCellInfo GmshCellTable[] = {
122: GmshCellEntry(15, VTX, 0, 0),
124: GmshCellEntry(1, SEG, 1, 1), GmshCellEntry(8, SEG, 1, 2), GmshCellEntry(26, SEG, 1, 3),
125: GmshCellEntry(27, SEG, 1, 4), GmshCellEntry(28, SEG, 1, 5), GmshCellEntry(62, SEG, 1, 6),
126: GmshCellEntry(63, SEG, 1, 7), GmshCellEntry(64, SEG, 1, 8), GmshCellEntry(65, SEG, 1, 9),
127: GmshCellEntry(66, SEG, 1, 10),
129: GmshCellEntry(2, TRI, 2, 1), GmshCellEntry(9, TRI, 2, 2), GmshCellEntry(21, TRI, 2, 3),
130: GmshCellEntry(23, TRI, 2, 4), GmshCellEntry(25, TRI, 2, 5), GmshCellEntry(42, TRI, 2, 6),
131: GmshCellEntry(43, TRI, 2, 7), GmshCellEntry(44, TRI, 2, 8), GmshCellEntry(45, TRI, 2, 9),
132: GmshCellEntry(46, TRI, 2, 10),
134: GmshCellEntry(3, QUA, 2, 1), GmshCellEntry(10, QUA, 2, 2), {16, GMSH_QUA, 2, 2, 4, 8, GmshLexOrder_QUA_2_Serendipity},
135: GmshCellEntry(36, QUA, 2, 3), GmshCellEntry(37, QUA, 2, 4), GmshCellEntry(38, QUA, 2, 5),
136: GmshCellEntry(47, QUA, 2, 6), GmshCellEntry(48, QUA, 2, 7), GmshCellEntry(49, QUA, 2, 8),
137: GmshCellEntry(50, QUA, 2, 9), GmshCellEntry(51, QUA, 2, 10),
139: GmshCellEntry(4, TET, 3, 1), GmshCellEntry(11, TET, 3, 2), GmshCellEntry(29, TET, 3, 3),
140: GmshCellEntry(30, TET, 3, 4), GmshCellEntry(31, TET, 3, 5), GmshCellEntry(71, TET, 3, 6),
141: GmshCellEntry(72, TET, 3, 7), GmshCellEntry(73, TET, 3, 8), GmshCellEntry(74, TET, 3, 9),
142: GmshCellEntry(75, TET, 3, 10),
144: GmshCellEntry(5, HEX, 3, 1), GmshCellEntry(12, HEX, 3, 2), {17, GMSH_HEX, 3, 2, 8, 20, GmshLexOrder_HEX_2_Serendipity},
145: GmshCellEntry(92, HEX, 3, 3), GmshCellEntry(93, HEX, 3, 4), GmshCellEntry(94, HEX, 3, 5),
146: GmshCellEntry(95, HEX, 3, 6), GmshCellEntry(96, HEX, 3, 7), GmshCellEntry(97, HEX, 3, 8),
147: GmshCellEntry(98, HEX, 3, 9), GmshCellEntry(-1, HEX, 3, 10),
149: GmshCellEntry(6, PRI, 3, 1), GmshCellEntry(13, PRI, 3, 2), GmshCellEntry(90, PRI, 3, 3),
150: GmshCellEntry(91, PRI, 3, 4), GmshCellEntry(106, PRI, 3, 5), GmshCellEntry(107, PRI, 3, 6),
151: GmshCellEntry(108, PRI, 3, 7), GmshCellEntry(109, PRI, 3, 8), GmshCellEntry(110, PRI, 3, 9),
152: GmshCellEntry(-1, PRI, 3, 10),
154: GmshCellEntry(7, PYR, 3, 1), GmshCellEntry(14, PYR, 3, 2), GmshCellEntry(118, PYR, 3, 3),
155: GmshCellEntry(119, PYR, 3, 4), GmshCellEntry(120, PYR, 3, 5), GmshCellEntry(121, PYR, 3, 6),
156: GmshCellEntry(122, PYR, 3, 7), GmshCellEntry(123, PYR, 3, 8), GmshCellEntry(124, PYR, 3, 9),
157: GmshCellEntry(-1, PYR, 3, 10)
159: #if 0
160: {20, GMSH_TRI, 2, 3, 3, 9, NULL},
161: {18, GMSH_PRI, 3, 2, 6, 15, NULL},
162: {19, GMSH_PYR, 3, 2, 5, 13, NULL},
163: #endif
164: };
166: static GmshCellInfo GmshCellMap[150];
168: static PetscErrorCode GmshCellInfoSetUp(void)
169: {
170: size_t i, n;
171: static PetscBool called = PETSC_FALSE;
173: PetscFunctionBegin;
174: if (called) PetscFunctionReturn(PETSC_SUCCESS);
175: called = PETSC_TRUE;
176: n = PETSC_STATIC_ARRAY_LENGTH(GmshCellMap);
177: for (i = 0; i < n; ++i) {
178: GmshCellMap[i].cellType = -1;
179: GmshCellMap[i].polytope = -1;
180: }
181: n = PETSC_STATIC_ARRAY_LENGTH(GmshCellTable);
182: for (i = 0; i < n; ++i) {
183: if (GmshCellTable[i].cellType <= 0) continue;
184: GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i];
185: }
186: PetscFunctionReturn(PETSC_SUCCESS);
187: }
189: #define GmshCellTypeCheck(ct) \
190: PetscMacroReturnStandard(const int _ct_ = (int)ct; PetscCheck(_ct_ >= 0 && _ct_ < (int)PETSC_STATIC_ARRAY_LENGTH(GmshCellMap), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); PetscCheck(GmshCellMap[_ct_].cellType == _ct_, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \
191: PetscCheck(GmshCellMap[_ct_].polytope != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_);)
193: typedef struct {
194: PetscViewer viewer;
195: int fileFormat;
196: int dataSize;
197: PetscBool binary;
198: PetscBool byteSwap;
199: size_t wlen;
200: void *wbuf;
201: size_t slen;
202: void *sbuf;
203: PetscInt *nbuf;
204: PetscInt nodeStart;
205: PetscInt nodeEnd;
206: PetscInt *nodeMap;
207: } GmshFile;
209: /*
210: Returns an array of count items each with a sizeof(eltsize)
211: */
212: static PetscErrorCode GmshBufferGet(GmshFile *gmsh, PetscCount count, size_t eltsize, void *buf)
213: {
214: size_t size = count * eltsize;
216: PetscFunctionBegin;
217: if (gmsh->wlen < size) {
218: PetscCall(PetscFree(gmsh->wbuf));
219: PetscCall(PetscMalloc(size, &gmsh->wbuf));
220: gmsh->wlen = size;
221: }
222: *(void **)buf = size ? gmsh->wbuf : NULL;
223: PetscFunctionReturn(PETSC_SUCCESS);
224: }
226: /*
227: Returns an array of count items each with the size determined by the GmshFile
228: */
229: static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, PetscCount count, void *buf)
230: {
231: size_t dataSize = (size_t)gmsh->dataSize;
232: size_t size = count * dataSize;
234: PetscFunctionBegin;
235: if (gmsh->slen < size) {
236: PetscCall(PetscFree(gmsh->sbuf));
237: PetscCall(PetscMalloc(size, &gmsh->sbuf));
238: gmsh->slen = size;
239: }
240: *(void **)buf = size ? gmsh->sbuf : NULL;
241: PetscFunctionReturn(PETSC_SUCCESS);
242: }
244: /*
245: Reads an array of count items each with the size determined by the PetscDataType
246: */
247: static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscCount count, PetscDataType dtype)
248: {
249: PetscInt icount;
251: PetscFunctionBegin;
252: PetscCall(PetscIntCast(count, &icount));
253: PetscCall(PetscViewerRead(gmsh->viewer, buf, icount, NULL, dtype));
254: if (gmsh->byteSwap) PetscCall(PetscByteSwap(buf, dtype, icount));
255: PetscFunctionReturn(PETSC_SUCCESS);
256: }
258: static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscCount count)
259: {
260: PetscInt icount;
262: PetscFunctionBegin;
263: PetscCall(PetscIntCast(count, &icount));
264: PetscCall(PetscViewerRead(gmsh->viewer, buf, icount, NULL, PETSC_STRING));
265: PetscFunctionReturn(PETSC_SUCCESS);
266: }
268: static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match)
269: {
270: PetscFunctionBegin;
271: PetscCall(PetscStrcmp(line, Section, match));
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN])
276: {
277: PetscBool match;
279: PetscFunctionBegin;
280: PetscCall(GmshMatch(gmsh, Section, line, &match));
281: PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s, not %s", Section, line);
282: PetscFunctionReturn(PETSC_SUCCESS);
283: }
285: static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN])
286: {
287: PetscBool match;
289: PetscFunctionBegin;
290: while (PETSC_TRUE) {
291: PetscCall(GmshReadString(gmsh, line, 1));
292: PetscCall(GmshMatch(gmsh, "$Comments", line, &match));
293: if (!match) break;
294: while (PETSC_TRUE) {
295: PetscCall(GmshReadString(gmsh, line, 1));
296: PetscCall(GmshMatch(gmsh, "$EndComments", line, &match));
297: if (match) break;
298: }
299: }
300: PetscFunctionReturn(PETSC_SUCCESS);
301: }
303: static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN])
304: {
305: PetscFunctionBegin;
306: PetscCall(GmshReadString(gmsh, line, 1));
307: PetscCall(GmshExpect(gmsh, EndSection, line));
308: PetscFunctionReturn(PETSC_SUCCESS);
309: }
311: /*
312: Read into buf[] count number of PetscInt integers (the file storage size may be different than PetscInt)
313: */
314: static PetscErrorCode GmshReadPetscInt(GmshFile *gmsh, PetscInt *buf, PetscCount count)
315: {
316: PetscCount i;
317: size_t dataSize = (size_t)gmsh->dataSize;
319: PetscFunctionBegin;
320: if (dataSize == sizeof(PetscInt)) PetscCall(GmshRead(gmsh, buf, count, PETSC_INT));
321: else if (dataSize == sizeof(int)) {
322: int *ibuf = NULL;
323: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
324: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_ENUM));
325: for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
326: } else if (dataSize == sizeof(long)) {
327: long *ibuf = NULL;
328: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
329: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_LONG));
330: for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
331: } else if (dataSize == sizeof(PetscInt64)) {
332: PetscInt64 *ibuf = NULL;
333: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
334: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_INT64));
335: for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i];
336: }
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: /*
341: Read into buf[] count number of PetscCount integers (the file storage size may be different than PetscCount)
342: */
343: static PetscErrorCode GmshReadPetscCount(GmshFile *gmsh, PetscCount *buf, PetscCount count)
344: {
345: PetscCount i;
346: size_t dataSize = (size_t)gmsh->dataSize;
348: PetscFunctionBegin;
349: if (dataSize == sizeof(PetscCount)) PetscCall(GmshRead(gmsh, buf, count, PETSC_COUNT));
350: else if (dataSize == sizeof(int)) {
351: int *ibuf = NULL;
352: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
353: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_ENUM));
354: for (i = 0; i < count; ++i) buf[i] = (PetscCount)ibuf[i];
355: } else if (dataSize == sizeof(long)) {
356: long *ibuf = NULL;
357: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
358: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_LONG));
359: for (i = 0; i < count; ++i) buf[i] = (PetscCount)ibuf[i];
360: } else if (dataSize == sizeof(PetscInt64)) {
361: PetscInt64 *ibuf = NULL;
362: PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf));
363: PetscCall(GmshRead(gmsh, ibuf, count, PETSC_INT64));
364: for (i = 0; i < count; ++i) buf[i] = (PetscCount)ibuf[i];
365: }
366: PetscFunctionReturn(PETSC_SUCCESS);
367: }
369: /*
370: Read into buf[] count number of PetscEnum integers
371: */
372: static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscCount count)
373: {
374: PetscFunctionBegin;
375: PetscCall(GmshRead(gmsh, buf, count, PETSC_ENUM));
376: PetscFunctionReturn(PETSC_SUCCESS);
377: }
379: /*
380: Read into buf[] count number of double
381: */
382: static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscCount count)
383: {
384: PetscFunctionBegin;
385: PetscCall(GmshRead(gmsh, buf, count, PETSC_DOUBLE));
386: PetscFunctionReturn(PETSC_SUCCESS);
387: }
389: #define GMSH_MAX_TAGS 16
391: typedef struct {
392: PetscInt id; /* Entity ID */
393: PetscInt dim; /* Dimension */
394: double bbox[6]; /* Bounding box */
395: PetscInt numTags; /* Size of tag array */
396: int tags[GMSH_MAX_TAGS]; /* Tag array */
397: } GmshEntity;
399: typedef struct {
400: GmshEntity *entity[4];
401: PetscHMapI entityMap[4];
402: } GmshEntities;
404: static PetscErrorCode GmshEntitiesCreate(PetscCount count[4], GmshEntities **entities)
405: {
406: PetscFunctionBegin;
407: PetscCall(PetscNew(entities));
408: for (PetscInt dim = 0; dim < 4; ++dim) {
409: PetscCall(PetscCalloc1(count[dim], &(*entities)->entity[dim]));
410: PetscCall(PetscHMapICreate(&(*entities)->entityMap[dim]));
411: }
412: PetscFunctionReturn(PETSC_SUCCESS);
413: }
415: static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities)
416: {
417: PetscInt dim;
419: PetscFunctionBegin;
420: if (!*entities) PetscFunctionReturn(PETSC_SUCCESS);
421: for (dim = 0; dim < 4; ++dim) {
422: PetscCall(PetscFree((*entities)->entity[dim]));
423: PetscCall(PetscHMapIDestroy(&(*entities)->entityMap[dim]));
424: }
425: PetscCall(PetscFree(*entities));
426: PetscFunctionReturn(PETSC_SUCCESS);
427: }
429: static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity **entity)
430: {
431: PetscFunctionBegin;
432: PetscCall(PetscHMapISet(entities->entityMap[dim], eid, index));
433: entities->entity[dim][index].dim = dim;
434: entities->entity[dim][index].id = eid;
435: if (entity) *entity = &entities->entity[dim][index];
436: PetscFunctionReturn(PETSC_SUCCESS);
437: }
439: static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity **entity)
440: {
441: PetscInt index;
443: PetscFunctionBegin;
444: PetscCall(PetscHMapIGet(entities->entityMap[dim], eid, &index));
445: *entity = &entities->entity[dim][index];
446: PetscFunctionReturn(PETSC_SUCCESS);
447: }
449: typedef struct {
450: PetscInt *id; /* Node IDs */
451: double *xyz; /* Coordinates */
452: PetscInt *tag; /* Physical tag */
453: } GmshNodes;
455: static PetscErrorCode GmshNodesCreate(PetscCount count, GmshNodes **nodes)
456: {
457: PetscFunctionBegin;
458: PetscCall(PetscNew(nodes));
459: PetscCall(PetscMalloc1(count * 1, &(*nodes)->id));
460: PetscCall(PetscMalloc1(count * 3, &(*nodes)->xyz));
461: PetscCall(PetscMalloc1(count * GMSH_MAX_TAGS, &(*nodes)->tag));
462: PetscFunctionReturn(PETSC_SUCCESS);
463: }
465: static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes)
466: {
467: PetscFunctionBegin;
468: if (!*nodes) PetscFunctionReturn(PETSC_SUCCESS);
469: PetscCall(PetscFree((*nodes)->id));
470: PetscCall(PetscFree((*nodes)->xyz));
471: PetscCall(PetscFree((*nodes)->tag));
472: PetscCall(PetscFree(*nodes));
473: PetscFunctionReturn(PETSC_SUCCESS);
474: }
476: typedef struct {
477: PetscInt id; /* Element ID */
478: PetscInt dim; /* Dimension */
479: PetscInt cellType; /* Cell type */
480: PetscInt numVerts; /* Size of vertex array */
481: PetscInt numNodes; /* Size of node array */
482: PetscInt *nodes; /* Vertex/Node array */
483: PetscInt numTags; /* Size of physical tag array */
484: int tags[GMSH_MAX_TAGS]; /* Physical tag array */
485: } GmshElement;
487: static PetscErrorCode GmshElementsCreate(PetscCount count, GmshElement **elements)
488: {
489: PetscFunctionBegin;
490: PetscCall(PetscCalloc1(count, elements));
491: PetscFunctionReturn(PETSC_SUCCESS);
492: }
494: static PetscErrorCode GmshElementsDestroy(GmshElement **elements)
495: {
496: PetscFunctionBegin;
497: if (!*elements) PetscFunctionReturn(PETSC_SUCCESS);
498: PetscCall(PetscFree(*elements));
499: PetscFunctionReturn(PETSC_SUCCESS);
500: }
502: typedef struct {
503: PetscInt dim;
504: PetscInt order;
505: GmshEntities *entities;
506: PetscInt numNodes;
507: GmshNodes *nodelist;
508: PetscInt numElems;
509: GmshElement *elements;
510: PetscInt numVerts;
511: PetscInt numCells;
512: PetscInt *periodMap;
513: PetscInt *vertexMap;
514: PetscSegBuffer segbuf;
515: PetscInt numRegions;
516: PetscInt *regionDims;
517: PetscInt *regionTags;
518: char **regionNames;
519: } GmshMesh;
521: static PetscErrorCode GmshMeshCreate(GmshMesh **mesh)
522: {
523: PetscFunctionBegin;
524: PetscCall(PetscNew(mesh));
525: PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf));
526: PetscFunctionReturn(PETSC_SUCCESS);
527: }
529: static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh)
530: {
531: PetscFunctionBegin;
532: if (!*mesh) PetscFunctionReturn(PETSC_SUCCESS);
533: PetscCall(GmshEntitiesDestroy(&(*mesh)->entities));
534: PetscCall(GmshNodesDestroy(&(*mesh)->nodelist));
535: PetscCall(GmshElementsDestroy(&(*mesh)->elements));
536: PetscCall(PetscFree((*mesh)->periodMap));
537: PetscCall(PetscFree((*mesh)->vertexMap));
538: PetscCall(PetscSegBufferDestroy(&(*mesh)->segbuf));
539: for (PetscInt r = 0; r < (*mesh)->numRegions; ++r) PetscCall(PetscFree((*mesh)->regionNames[r]));
540: PetscCall(PetscFree3((*mesh)->regionDims, (*mesh)->regionTags, (*mesh)->regionNames));
541: PetscCall(PetscFree(*mesh));
542: PetscFunctionReturn(PETSC_SUCCESS);
543: }
545: static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh)
546: {
547: PetscViewer viewer = gmsh->viewer;
548: PetscBool byteSwap = gmsh->byteSwap;
549: char line[PETSC_MAX_PATH_LEN];
550: int n, t, num, nid, snum;
551: GmshNodes *nodes;
553: PetscFunctionBegin;
554: PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
555: snum = sscanf(line, "%d", &num);
556: PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
557: PetscCall(GmshNodesCreate(num, &nodes));
558: mesh->numNodes = num;
559: mesh->nodelist = nodes;
560: for (n = 0; n < num; ++n) {
561: double *xyz = nodes->xyz + n * 3;
562: PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM));
563: PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE));
564: if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
565: if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
566: nodes->id[n] = nid;
567: for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1;
568: }
569: PetscFunctionReturn(PETSC_SUCCESS);
570: }
572: /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the
573: file contents multiple times to figure out the true number of cells and facets
574: in the given mesh. To make this more efficient we read the file contents only
575: once and store them in memory, while determining the true number of cells. */
576: static PetscErrorCode GmshReadElements_v22(GmshFile *gmsh, GmshMesh *mesh)
577: {
578: PetscViewer viewer = gmsh->viewer;
579: PetscBool binary = gmsh->binary;
580: PetscBool byteSwap = gmsh->byteSwap;
581: char line[PETSC_MAX_PATH_LEN];
582: int i, c, p, num, ibuf[1 + 4 + 1000], snum;
583: int cellType, numElem, numVerts, numNodes, numTags;
584: GmshElement *elements;
585: PetscInt *nodeMap = gmsh->nodeMap;
587: PetscFunctionBegin;
588: PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
589: snum = sscanf(line, "%d", &num);
590: PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
591: PetscCall(GmshElementsCreate(num, &elements));
592: mesh->numElems = num;
593: mesh->elements = elements;
594: for (c = 0; c < num;) {
595: PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM));
596: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3));
598: cellType = binary ? ibuf[0] : ibuf[1];
599: numElem = binary ? ibuf[1] : 1;
600: numTags = ibuf[2];
602: PetscCall(GmshCellTypeCheck(cellType));
603: numVerts = GmshCellMap[cellType].numVerts;
604: numNodes = GmshCellMap[cellType].numNodes;
606: for (i = 0; i < numElem; ++i, ++c) {
607: GmshElement *element = elements + c;
608: const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off;
609: const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1;
610: PetscCall(PetscViewerRead(viewer, ibuf + off, nint, NULL, PETSC_ENUM));
611: if (byteSwap) PetscCall(PetscByteSwap(ibuf + off, PETSC_ENUM, nint));
612: element->id = id[0];
613: element->dim = GmshCellMap[cellType].dim;
614: element->cellType = cellType;
615: element->numVerts = numVerts;
616: element->numNodes = numNodes;
617: element->numTags = PetscMin(numTags, GMSH_MAX_TAGS);
618: PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes));
619: for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
620: for (p = 0; p < element->numTags; p++) element->tags[p] = tags[p];
621: }
622: }
623: PetscFunctionReturn(PETSC_SUCCESS);
624: }
626: /*
627: $Entities
628: numPoints(unsigned long) numCurves(unsigned long)
629: numSurfaces(unsigned long) numVolumes(unsigned long)
630: // points
631: tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
632: boxMaxX(double) boxMaxY(double) boxMaxZ(double)
633: numPhysicals(unsigned long) phyisicalTag[...](int)
634: ...
635: // curves
636: tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
637: boxMaxX(double) boxMaxY(double) boxMaxZ(double)
638: numPhysicals(unsigned long) physicalTag[...](int)
639: numBREPVert(unsigned long) tagBREPVert[...](int)
640: ...
641: // surfaces
642: tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
643: boxMaxX(double) boxMaxY(double) boxMaxZ(double)
644: numPhysicals(unsigned long) physicalTag[...](int)
645: numBREPCurve(unsigned long) tagBREPCurve[...](int)
646: ...
647: // volumes
648: tag(int) boxMinX(double) boxMinY(double) boxMinZ(double)
649: boxMaxX(double) boxMaxY(double) boxMaxZ(double)
650: numPhysicals(unsigned long) physicalTag[...](int)
651: numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int)
652: ...
653: $EndEntities
654: */
655: static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh)
656: {
657: PetscViewer viewer = gmsh->viewer;
658: PetscBool byteSwap = gmsh->byteSwap;
659: long lnum, lbuf[4];
660: int dim, eid, numTags, *ibuf, t;
661: PetscCount index, count[4];
662: PetscInt i, num;
663: GmshEntity *entity = NULL;
665: PetscFunctionBegin;
666: PetscCall(PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG));
667: if (byteSwap) PetscCall(PetscByteSwap(lbuf, PETSC_LONG, 4));
668: for (i = 0; i < 4; ++i) count[i] = (PetscCount)lbuf[i];
669: PetscCall(GmshEntitiesCreate(count, &mesh->entities));
670: for (dim = 0; dim < 4; ++dim) {
671: for (index = 0; index < count[dim]; ++index) {
672: PetscCall(PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM));
673: if (byteSwap) PetscCall(PetscByteSwap(&eid, PETSC_ENUM, 1));
674: PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity));
675: PetscCall(PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE));
676: if (byteSwap) PetscCall(PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6));
677: PetscCall(PetscViewerRead(viewer, &lnum, 1, NULL, PETSC_LONG));
678: if (byteSwap) PetscCall(PetscByteSwap(&lnum, PETSC_LONG, 1));
679: PetscCall(PetscIntCast(lnum, &num));
680: PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf));
681: PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM));
682: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num));
683: entity->numTags = numTags = (int)PetscMin(num, GMSH_MAX_TAGS);
684: for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t];
685: if (dim == 0) continue;
686: PetscCall(PetscViewerRead(viewer, &lnum, 1, NULL, PETSC_LONG));
687: if (byteSwap) PetscCall(PetscByteSwap(&lnum, PETSC_LONG, 1));
688: PetscCall(GmshBufferGet(gmsh, lnum, sizeof(int), &ibuf));
689: PetscCall(PetscIntCast(lnum, &num));
690: PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM));
691: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num));
692: }
693: }
694: PetscFunctionReturn(PETSC_SUCCESS);
695: }
697: /*
698: $Nodes
699: numEntityBlocks(unsigned long) numNodes(unsigned long)
700: tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long)
701: tag(int) x(double) y(double) z(double)
702: ...
703: ...
704: $EndNodes
705: */
706: static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh)
707: {
708: PetscViewer viewer = gmsh->viewer;
709: PetscBool byteSwap = gmsh->byteSwap;
710: long block, node, n, t, numEntityBlocks, numTotalNodes, numNodes;
711: int info[3], nid;
712: GmshNodes *nodes;
714: PetscFunctionBegin;
715: PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG));
716: if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1));
717: PetscCall(PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG));
718: if (byteSwap) PetscCall(PetscByteSwap(&numTotalNodes, PETSC_LONG, 1));
719: PetscCall(GmshNodesCreate(numTotalNodes, &nodes));
720: PetscCall(PetscIntCast(numTotalNodes, &mesh->numNodes));
721: mesh->nodelist = nodes;
722: for (n = 0, block = 0; block < numEntityBlocks; ++block) {
723: PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM));
724: PetscCall(PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG));
725: if (byteSwap) PetscCall(PetscByteSwap(&numNodes, PETSC_LONG, 1));
726: if (gmsh->binary) {
727: size_t nbytes = sizeof(int) + 3 * sizeof(double);
728: char *cbuf = NULL; /* dummy value to prevent warning from compiler about possible uninitialized value */
729: PetscInt icnt;
731: PetscCall(GmshBufferGet(gmsh, numNodes, nbytes, &cbuf));
732: PetscCall(PetscIntCast(numNodes * nbytes, &icnt));
733: PetscCall(PetscViewerRead(viewer, cbuf, icnt, NULL, PETSC_CHAR));
734: for (node = 0; node < numNodes; ++node, ++n) {
735: char *cnid = cbuf + node * nbytes, *cxyz = cnid + sizeof(int);
736: double *xyz = nodes->xyz + n * 3;
738: if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cnid, PETSC_ENUM, 1));
739: if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cxyz, PETSC_DOUBLE, 3));
740: PetscCall(PetscMemcpy(&nid, cnid, sizeof(int)));
741: PetscCall(PetscMemcpy(xyz, cxyz, 3 * sizeof(double)));
742: if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
743: if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
744: nodes->id[n] = nid;
745: for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1;
746: }
747: } else {
748: for (node = 0; node < numNodes; ++node, ++n) {
749: double *xyz = nodes->xyz + n * 3;
751: PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM));
752: PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE));
753: if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1));
754: if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3));
755: nodes->id[n] = nid;
756: for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1;
757: }
758: }
759: }
760: PetscFunctionReturn(PETSC_SUCCESS);
761: }
763: /*
764: $Elements
765: numEntityBlocks(unsigned long) numElements(unsigned long)
766: tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long)
767: tag(int) numVert[...](int)
768: ...
769: ...
770: $EndElements
771: */
772: static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh)
773: {
774: PetscViewer viewer = gmsh->viewer;
775: PetscBool byteSwap = gmsh->byteSwap;
776: long c, block, numEntityBlocks, numTotalElements, elem, numElements;
777: int p, info[3], *ibuf = NULL;
778: int eid, dim, cellType, numVerts, numNodes, numTags;
779: GmshEntity *entity = NULL;
780: GmshElement *elements;
781: PetscInt *nodeMap = gmsh->nodeMap, icnt;
783: PetscFunctionBegin;
784: PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG));
785: if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1));
786: PetscCall(PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG));
787: if (byteSwap) PetscCall(PetscByteSwap(&numTotalElements, PETSC_LONG, 1));
788: PetscCall(GmshElementsCreate(numTotalElements, &elements));
789: PetscCall(PetscIntCast(numTotalElements, &mesh->numElems));
790: mesh->elements = elements;
791: for (c = 0, block = 0; block < numEntityBlocks; ++block) {
792: PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM));
793: if (byteSwap) PetscCall(PetscByteSwap(info, PETSC_ENUM, 3));
794: eid = info[0];
795: dim = info[1];
796: cellType = info[2];
797: PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
798: PetscCall(GmshCellTypeCheck(cellType));
799: numVerts = GmshCellMap[cellType].numVerts;
800: numNodes = GmshCellMap[cellType].numNodes;
801: numTags = (int)entity->numTags;
802: PetscCall(PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG));
803: if (byteSwap) PetscCall(PetscByteSwap(&numElements, PETSC_LONG, 1));
804: PetscCall(GmshBufferGet(gmsh, (1 + numNodes) * numElements, sizeof(int), &ibuf));
805: PetscCall(PetscIntCast((1 + numNodes) * numElements, &icnt));
806: PetscCall(PetscViewerRead(viewer, ibuf, icnt, NULL, PETSC_ENUM));
807: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, (1 + numNodes) * numElements));
808: for (elem = 0; elem < numElements; ++elem, ++c) {
809: GmshElement *element = elements + c;
810: const int *id = ibuf + elem * (1 + numNodes), *nodes = id + 1;
811: element->id = id[0];
812: element->dim = dim;
813: element->cellType = cellType;
814: element->numVerts = numVerts;
815: element->numNodes = numNodes;
816: element->numTags = numTags;
817: PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes));
818: for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
819: for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p];
820: }
821: }
822: PetscFunctionReturn(PETSC_SUCCESS);
823: }
825: static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[])
826: {
827: PetscViewer viewer = gmsh->viewer;
828: int fileFormat = gmsh->fileFormat;
829: PetscBool binary = gmsh->binary;
830: PetscBool byteSwap = gmsh->byteSwap;
831: int numPeriodic, snum, i;
832: char line[PETSC_MAX_PATH_LEN];
833: PetscInt *nodeMap = gmsh->nodeMap;
835: PetscFunctionBegin;
836: if (fileFormat == 22 || !binary) {
837: PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
838: snum = sscanf(line, "%d", &numPeriodic);
839: PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
840: } else {
841: PetscCall(PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM));
842: if (byteSwap) PetscCall(PetscByteSwap(&numPeriodic, PETSC_ENUM, 1));
843: }
844: for (i = 0; i < numPeriodic; i++) {
845: int ibuf[3], correspondingDim = -1, correspondingTag = -1, primaryTag = -1, correspondingNode, primaryNode;
846: long j, nNodes;
847: double affine[16];
849: if (fileFormat == 22 || !binary) {
850: PetscCall(PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING));
851: snum = sscanf(line, "%d %d %d", &correspondingDim, &correspondingTag, &primaryTag);
852: PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
853: } else {
854: PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM));
855: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3));
856: correspondingDim = ibuf[0];
857: correspondingTag = ibuf[1];
858: primaryTag = ibuf[2];
859: }
860: (void)correspondingDim;
861: (void)correspondingTag;
862: (void)primaryTag; /* unused */
864: if (fileFormat == 22 || !binary) {
865: PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
866: snum = sscanf(line, "%ld", &nNodes);
867: if (snum != 1) { /* discard transformation and try again */
868: PetscCall(PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING));
869: PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
870: snum = sscanf(line, "%ld", &nNodes);
871: PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
872: }
873: } else {
874: PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG));
875: if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1));
876: if (nNodes == -1) { /* discard transformation and try again */
877: PetscCall(PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE));
878: PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG));
879: if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1));
880: }
881: }
883: for (j = 0; j < nNodes; j++) {
884: if (fileFormat == 22 || !binary) {
885: PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING));
886: snum = sscanf(line, "%d %d", &correspondingNode, &primaryNode);
887: PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
888: } else {
889: PetscCall(PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM));
890: if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 2));
891: correspondingNode = ibuf[0];
892: primaryNode = ibuf[1];
893: }
894: correspondingNode = (int)nodeMap[correspondingNode];
895: primaryNode = (int)nodeMap[primaryNode];
896: periodicMap[correspondingNode] = primaryNode;
897: }
898: }
899: PetscFunctionReturn(PETSC_SUCCESS);
900: }
902: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
903: $Entities
904: numPoints(size_t) numCurves(size_t)
905: numSurfaces(size_t) numVolumes(size_t)
906: pointTag(int) X(double) Y(double) Z(double)
907: numPhysicalTags(size_t) physicalTag(int) ...
908: ...
909: curveTag(int) minX(double) minY(double) minZ(double)
910: maxX(double) maxY(double) maxZ(double)
911: numPhysicalTags(size_t) physicalTag(int) ...
912: numBoundingPoints(size_t) pointTag(int) ...
913: ...
914: surfaceTag(int) minX(double) minY(double) minZ(double)
915: maxX(double) maxY(double) maxZ(double)
916: numPhysicalTags(size_t) physicalTag(int) ...
917: numBoundingCurves(size_t) curveTag(int) ...
918: ...
919: volumeTag(int) minX(double) minY(double) minZ(double)
920: maxX(double) maxY(double) maxZ(double)
921: numPhysicalTags(size_t) physicalTag(int) ...
922: numBoundngSurfaces(size_t) surfaceTag(int) ...
923: ...
924: $EndEntities
925: */
926: static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh)
927: {
928: PetscCount count[4], index, numTags;
929: int dim, eid, *tags = NULL;
930: GmshEntity *entity = NULL;
932: PetscFunctionBegin;
933: PetscCall(GmshReadPetscCount(gmsh, count, 4));
934: PetscCall(GmshEntitiesCreate(count, &mesh->entities));
935: for (dim = 0; dim < 4; ++dim) {
936: for (index = 0; index < count[dim]; ++index) {
937: PetscCall(GmshReadInt(gmsh, &eid, 1));
938: PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity));
939: PetscCall(GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6));
940: PetscCall(GmshReadPetscCount(gmsh, &numTags, 1));
941: PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags));
942: PetscCall(GmshReadInt(gmsh, tags, numTags));
943: PetscCheck(numTags <= GMSH_MAX_TAGS, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PETSc currently supports up to %" PetscInt_FMT " tags per entity, not %" PetscCount_FMT, (PetscInt)GMSH_MAX_TAGS, numTags);
944: PetscCall(PetscIntCast(numTags, &entity->numTags));
945: for (PetscInt i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i];
946: if (dim == 0) continue;
947: PetscCall(GmshReadPetscCount(gmsh, &numTags, 1));
948: PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags));
949: PetscCall(GmshReadInt(gmsh, tags, numTags));
950: /* Currently, we do not save the ids for the bounding entities */
951: }
952: }
953: PetscFunctionReturn(PETSC_SUCCESS);
954: }
956: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
957: $Nodes
958: numEntityBlocks(size_t) numNodes(size_t)
959: minNodeTag(size_t) maxNodeTag(size_t)
960: entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t)
961: nodeTag(size_t)
962: ...
963: x(double) y(double) z(double)
964: < u(double; if parametric and entityDim = 1 or entityDim = 2) >
965: < v(double; if parametric and entityDim = 2) >
966: ...
967: ...
968: $EndNodes
969: */
970: static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh)
971: {
972: int info[3], dim, eid, parametric;
973: PetscCount sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0;
974: PetscInt numTags;
975: GmshEntity *entity = NULL;
976: GmshNodes *nodes;
978: PetscFunctionBegin;
979: PetscCall(GmshReadPetscCount(gmsh, sizes, 4));
980: numEntityBlocks = sizes[0];
981: numNodes = sizes[1];
982: PetscCall(GmshNodesCreate(numNodes, &nodes));
983: PetscCall(PetscIntCast(numNodes, &mesh->numNodes));
984: mesh->nodelist = nodes;
985: if (numEntityBlocks && !mesh->entities) PetscCall(PetscInfo(NULL, "File specifies %" PetscCount_FMT " entity blocks, but was missing the $Entities section\n", numEntityBlocks));
986: for (PetscCount block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) {
987: PetscCall(GmshReadInt(gmsh, info, 3));
988: dim = info[0];
989: eid = info[1];
990: parametric = info[2];
991: if (mesh->entities) PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
992: numTags = entity ? entity->numTags : 0;
993: PetscCheck(!parametric, PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported");
994: PetscCall(GmshReadPetscCount(gmsh, &numNodesBlock, 1));
995: PetscCall(GmshReadPetscInt(gmsh, nodes->id + node, numNodesBlock));
996: PetscCall(GmshReadDouble(gmsh, nodes->xyz + node * 3, numNodesBlock * 3));
997: for (PetscCount n = 0; n < numNodesBlock; ++n) {
998: PetscInt *tags = &nodes->tag[node * GMSH_MAX_TAGS];
1000: for (PetscInt t = 0; t < numTags; ++t) tags[n * GMSH_MAX_TAGS + t] = entity->tags[t];
1001: for (PetscInt t = numTags; t < GMSH_MAX_TAGS; ++t) tags[n * GMSH_MAX_TAGS + t] = -1;
1002: }
1003: }
1004: PetscCall(PetscIntCast(sizes[2], &gmsh->nodeStart));
1005: PetscCall(PetscIntCast(sizes[3] + 1, &gmsh->nodeEnd));
1006: PetscFunctionReturn(PETSC_SUCCESS);
1007: }
1009: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1010: $Elements
1011: numEntityBlocks(size_t) numElements(size_t)
1012: minElementTag(size_t) maxElementTag(size_t)
1013: entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t)
1014: elementTag(size_t) nodeTag(size_t) ...
1015: ...
1016: ...
1017: $EndElements
1018: */
1019: static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh)
1020: {
1021: int info[3], eid, dim, cellType;
1022: PetscCount sizes[4], numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p;
1023: GmshEntity *entity = NULL;
1024: GmshElement *elements;
1025: PetscInt *nodeMap = gmsh->nodeMap, *ibuf = NULL;
1027: PetscFunctionBegin;
1028: PetscCall(GmshReadPetscCount(gmsh, sizes, 4));
1029: numEntityBlocks = sizes[0];
1030: numElements = sizes[1];
1031: PetscCall(GmshElementsCreate(numElements, &elements));
1032: PetscCall(PetscIntCast(numElements, &mesh->numElems));
1033: mesh->elements = elements;
1034: if (numEntityBlocks && !mesh->entities) PetscCall(PetscInfo(NULL, "File specifies %" PetscCount_FMT " entity blocks, but was missing the $Entities section\n", numEntityBlocks));
1035: for (c = 0, block = 0; block < numEntityBlocks; ++block) {
1036: PetscCall(GmshReadInt(gmsh, info, 3));
1037: dim = info[0];
1038: eid = info[1];
1039: cellType = info[2];
1040: if (mesh->entities) PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity));
1041: PetscCall(GmshCellTypeCheck(cellType));
1042: numVerts = GmshCellMap[cellType].numVerts;
1043: numNodes = GmshCellMap[cellType].numNodes;
1044: numTags = entity ? entity->numTags : 0;
1045: PetscCall(GmshReadPetscCount(gmsh, &numBlockElements, 1));
1046: PetscCall(GmshBufferGet(gmsh, (1 + numNodes) * numBlockElements, sizeof(PetscInt), &ibuf));
1047: PetscCall(GmshReadPetscInt(gmsh, ibuf, (1 + numNodes) * numBlockElements));
1048: for (elem = 0; elem < numBlockElements; ++elem, ++c) {
1049: GmshElement *element = elements + c;
1050: const PetscInt *id = ibuf + elem * (1 + numNodes), *nodes = id + 1;
1052: element->id = id[0];
1053: element->dim = dim;
1054: element->cellType = cellType;
1055: PetscCall(PetscIntCast(numVerts, &element->numVerts));
1056: PetscCall(PetscIntCast(numNodes, &element->numNodes));
1057: PetscCall(PetscIntCast(numTags, &element->numTags));
1058: PetscCall(PetscSegBufferGet(mesh->segbuf, element->numNodes, &element->nodes));
1059: for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]];
1060: for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p];
1061: }
1062: }
1063: PetscFunctionReturn(PETSC_SUCCESS);
1064: }
1066: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1067: $Periodic
1068: numPeriodicLinks(size_t)
1069: entityDim(int) entityTag(int) entityTagPrimary(int)
1070: numAffine(size_t) value(double) ...
1071: numCorrespondingNodes(size_t)
1072: nodeTag(size_t) nodeTagPrimary(size_t)
1073: ...
1074: ...
1075: $EndPeriodic
1076: */
1077: static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[])
1078: {
1079: int info[3];
1080: double dbuf[16];
1081: PetscCount numPeriodicLinks, numAffine, numCorrespondingNodes;
1082: PetscInt *nodeMap = gmsh->nodeMap, *nodeTags = NULL;
1084: PetscFunctionBegin;
1085: PetscCall(GmshReadPetscCount(gmsh, &numPeriodicLinks, 1));
1086: for (PetscCount link = 0; link < numPeriodicLinks; ++link) {
1087: PetscCall(GmshReadInt(gmsh, info, 3));
1088: PetscCall(GmshReadPetscCount(gmsh, &numAffine, 1));
1089: PetscCall(GmshReadDouble(gmsh, dbuf, numAffine));
1090: PetscCall(GmshReadPetscCount(gmsh, &numCorrespondingNodes, 1));
1091: PetscCall(GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags));
1092: PetscCall(GmshReadPetscInt(gmsh, nodeTags, numCorrespondingNodes * 2));
1093: for (PetscCount node = 0; node < numCorrespondingNodes; ++node) {
1094: PetscInt correspondingNode = nodeMap[nodeTags[node * 2 + 0]];
1095: PetscInt primaryNode = nodeMap[nodeTags[node * 2 + 1]];
1097: periodicMap[correspondingNode] = primaryNode;
1098: }
1099: }
1100: PetscFunctionReturn(PETSC_SUCCESS);
1101: }
1103: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1104: $MeshFormat // same as MSH version 2
1105: version(ASCII double; currently 4.1)
1106: file-type(ASCII int; 0 for ASCII mode, 1 for binary mode)
1107: data-size(ASCII int; sizeof(size_t))
1108: < int with value one; only in binary mode, to detect endianness >
1109: $EndMeshFormat
1110: */
1111: static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh)
1112: {
1113: char line[PETSC_MAX_PATH_LEN];
1114: int snum, fileType, fileFormat, dataSize, checkEndian;
1115: float version;
1117: PetscFunctionBegin;
1118: PetscCall(GmshReadString(gmsh, line, 3));
1119: snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize);
1120: fileFormat = (int)roundf(version * 10);
1121: PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
1122: PetscCheck(fileFormat >= 22, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
1123: PetscCheck((int)version != 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
1124: PetscCheck(fileFormat <= 41, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
1125: PetscCheck(!gmsh->binary || fileType, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII");
1126: PetscCheck(gmsh->binary || !fileType, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary");
1127: PetscCheck(fileFormat > 40 || dataSize == sizeof(double), PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize);
1128: PetscCheck(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);
1129: gmsh->fileFormat = fileFormat;
1130: gmsh->dataSize = dataSize;
1131: gmsh->byteSwap = PETSC_FALSE;
1132: if (gmsh->binary) {
1133: PetscCall(GmshReadInt(gmsh, &checkEndian, 1));
1134: if (checkEndian != 1) {
1135: PetscCall(PetscByteSwap(&checkEndian, PETSC_ENUM, 1));
1136: PetscCheck(checkEndian == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line);
1137: gmsh->byteSwap = PETSC_TRUE;
1138: }
1139: }
1140: PetscFunctionReturn(PETSC_SUCCESS);
1141: }
1143: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1144: Neper: https://neper.info/ adds this section
1146: $MeshVersion
1147: <major>.<minor>,<patch>
1148: $EndMeshVersion
1149: */
1150: static PetscErrorCode GmshReadMeshVersion(GmshFile *gmsh)
1151: {
1152: char line[PETSC_MAX_PATH_LEN];
1153: int snum, major, minor, patch;
1155: PetscFunctionBegin;
1156: PetscCall(GmshReadString(gmsh, line, 1));
1157: snum = sscanf(line, "%d.%d.%d", &major, &minor, &patch);
1158: PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
1159: PetscFunctionReturn(PETSC_SUCCESS);
1160: }
1162: /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format
1163: Neper: https://neper.info/ adds this section
1165: $Domain
1166: <shape>
1167: $EndDomain
1168: */
1169: static PetscErrorCode GmshReadMeshDomain(GmshFile *gmsh)
1170: {
1171: char line[PETSC_MAX_PATH_LEN];
1173: PetscFunctionBegin;
1174: PetscCall(GmshReadString(gmsh, line, 1));
1175: PetscFunctionReturn(PETSC_SUCCESS);
1176: }
1178: /*
1179: PhysicalNames
1180: numPhysicalNames(ASCII int)
1181: dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max)
1182: ...
1183: $EndPhysicalNames
1184: */
1185: static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh, GmshMesh *mesh)
1186: {
1187: char line[PETSC_MAX_PATH_LEN], name[128 + 2], *p = NULL, *q = NULL, *r = NULL;
1188: int snum, region, dim, tag;
1190: PetscFunctionBegin;
1191: PetscCall(GmshReadString(gmsh, line, 1));
1192: snum = sscanf(line, "%d", ®ion);
1193: mesh->numRegions = region;
1194: PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1195: PetscCall(PetscMalloc3(mesh->numRegions, &mesh->regionDims, mesh->numRegions, &mesh->regionTags, mesh->numRegions, &mesh->regionNames));
1196: for (region = 0; region < mesh->numRegions; ++region) {
1197: PetscCall(GmshReadString(gmsh, line, 2));
1198: snum = sscanf(line, "%d %d", &dim, &tag);
1199: PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1200: PetscCall(GmshReadString(gmsh, line, -(PetscInt)sizeof(line)));
1201: PetscCall(PetscStrchr(line, '"', &p));
1202: PetscCheck(p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1203: PetscCall(PetscStrrchr(line, '"', &q));
1204: PetscCheck(q != p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file");
1205: PetscCall(PetscStrrchr(line, ':', &r));
1206: if (p != r) q = r;
1207: PetscCall(PetscStrncpy(name, p + 1, (size_t)(q - p - 1)));
1208: mesh->regionDims[region] = dim;
1209: mesh->regionTags[region] = tag;
1210: PetscCall(PetscStrallocpy(name, &mesh->regionNames[region]));
1211: }
1212: PetscFunctionReturn(PETSC_SUCCESS);
1213: }
1215: static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh)
1216: {
1217: PetscFunctionBegin;
1218: switch (gmsh->fileFormat) {
1219: case 41:
1220: PetscCall(GmshReadEntities_v41(gmsh, mesh));
1221: break;
1222: default:
1223: PetscCall(GmshReadEntities_v40(gmsh, mesh));
1224: break;
1225: }
1226: PetscFunctionReturn(PETSC_SUCCESS);
1227: }
1229: static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh)
1230: {
1231: PetscFunctionBegin;
1232: switch (gmsh->fileFormat) {
1233: case 41:
1234: PetscCall(GmshReadNodes_v41(gmsh, mesh));
1235: break;
1236: case 40:
1237: PetscCall(GmshReadNodes_v40(gmsh, mesh));
1238: break;
1239: default:
1240: PetscCall(GmshReadNodes_v22(gmsh, mesh));
1241: break;
1242: }
1244: { /* Gmsh v2.2/v4.0 does not provide min/max node tags */
1245: if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) {
1246: PetscInt tagMin = PETSC_INT_MAX, tagMax = PETSC_INT_MIN, n;
1247: GmshNodes *nodes = mesh->nodelist;
1248: for (n = 0; n < mesh->numNodes; ++n) {
1249: const PetscInt tag = nodes->id[n];
1250: tagMin = PetscMin(tag, tagMin);
1251: tagMax = PetscMax(tag, tagMax);
1252: }
1253: gmsh->nodeStart = tagMin;
1254: gmsh->nodeEnd = tagMax + 1;
1255: }
1256: }
1258: { /* Support for sparse node tags */
1259: PetscInt n;
1260: GmshNodes *nodes = mesh->nodelist;
1261: PetscCall(PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf));
1262: for (PetscInt t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_INT_MIN;
1263: gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart;
1264: for (n = 0; n < mesh->numNodes; ++n) {
1265: const PetscInt tag = nodes->id[n];
1266: PetscCheck(gmsh->nodeMap[tag] < 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %" PetscInt_FMT, tag);
1267: gmsh->nodeMap[tag] = n;
1268: }
1269: }
1270: PetscFunctionReturn(PETSC_SUCCESS);
1271: }
1273: static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh)
1274: {
1275: PetscFunctionBegin;
1276: switch (gmsh->fileFormat) {
1277: case 41:
1278: PetscCall(GmshReadElements_v41(gmsh, mesh));
1279: break;
1280: case 40:
1281: PetscCall(GmshReadElements_v40(gmsh, mesh));
1282: break;
1283: default:
1284: PetscCall(GmshReadElements_v22(gmsh, mesh));
1285: break;
1286: }
1288: { /* Reorder elements by codimension and polytope type */
1289: PetscInt ne = mesh->numElems;
1290: GmshElement *elements = mesh->elements;
1291: PetscInt keymap[GMSH_NUM_POLYTOPES], nk = 0;
1292: PetscInt offset[GMSH_NUM_POLYTOPES + 1], e, k;
1294: for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_INT_MIN;
1295: PetscCall(PetscMemzero(offset, sizeof(offset)));
1297: keymap[GMSH_TET] = nk++;
1298: keymap[GMSH_HEX] = nk++;
1299: keymap[GMSH_PRI] = nk++;
1300: keymap[GMSH_PYR] = nk++;
1301: keymap[GMSH_TRI] = nk++;
1302: keymap[GMSH_QUA] = nk++;
1303: keymap[GMSH_SEG] = nk++;
1304: keymap[GMSH_VTX] = nk++;
1306: PetscCall(GmshElementsCreate(mesh->numElems, &mesh->elements));
1307: #define key(eid) keymap[GmshCellMap[elements[eid].cellType].polytope]
1308: for (e = 0; e < ne; ++e) offset[1 + key(e)]++;
1309: for (k = 1; k < nk; ++k) offset[k] += offset[k - 1];
1310: for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e];
1311: #undef key
1312: PetscCall(GmshElementsDestroy(&elements));
1313: }
1315: { /* Mesh dimension and order */
1316: GmshElement *elem = mesh->numElems ? mesh->elements : NULL;
1317: mesh->dim = elem ? GmshCellMap[elem->cellType].dim : 0;
1318: mesh->order = elem ? GmshCellMap[elem->cellType].order : 0;
1319: }
1321: {
1322: PetscBT vtx;
1323: PetscInt dim = mesh->dim, e, n, v;
1325: PetscCall(PetscBTCreate(mesh->numNodes, &vtx));
1327: /* Compute number of cells and set of vertices */
1328: mesh->numCells = 0;
1329: for (e = 0; e < mesh->numElems; ++e) {
1330: GmshElement *elem = mesh->elements + e;
1331: if (elem->dim == dim && dim > 0) mesh->numCells++;
1332: for (v = 0; v < elem->numVerts; v++) PetscCall(PetscBTSet(vtx, elem->nodes[v]));
1333: }
1335: /* Compute numbering for vertices */
1336: mesh->numVerts = 0;
1337: PetscCall(PetscMalloc1(mesh->numNodes, &mesh->vertexMap));
1338: for (n = 0; n < mesh->numNodes; ++n) mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_INT_MIN;
1340: PetscCall(PetscBTDestroy(&vtx));
1341: }
1342: PetscFunctionReturn(PETSC_SUCCESS);
1343: }
1345: static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh)
1346: {
1347: PetscInt n;
1349: PetscFunctionBegin;
1350: PetscCall(PetscMalloc1(mesh->numNodes, &mesh->periodMap));
1351: for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n;
1352: switch (gmsh->fileFormat) {
1353: case 41:
1354: PetscCall(GmshReadPeriodic_v41(gmsh, mesh->periodMap));
1355: break;
1356: default:
1357: PetscCall(GmshReadPeriodic_v40(gmsh, mesh->periodMap));
1358: break;
1359: }
1361: /* Find canonical primary nodes */
1362: for (n = 0; n < mesh->numNodes; ++n)
1363: while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]]) mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]];
1365: /* Renumber vertices (filter out correspondings) */
1366: mesh->numVerts = 0;
1367: for (n = 0; n < mesh->numNodes; ++n)
1368: if (mesh->vertexMap[n] >= 0) /* is vertex */
1369: if (mesh->periodMap[n] == n) /* is primary */
1370: mesh->vertexMap[n] = mesh->numVerts++;
1371: for (n = 0; n < mesh->numNodes; ++n)
1372: if (mesh->vertexMap[n] >= 0) /* is vertex */
1373: if (mesh->periodMap[n] != n) /* is corresponding */
1374: mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]];
1375: PetscFunctionReturn(PETSC_SUCCESS);
1376: }
1378: #define DM_POLYTOPE_VERTEX DM_POLYTOPE_POINT
1379: static const DMPolytopeType DMPolytopeMap[] = {
1380: /* GMSH_VTX */ DM_POLYTOPE_VERTEX,
1381: /* GMSH_SEG */ DM_POLYTOPE_SEGMENT,
1382: /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE,
1383: /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL,
1384: /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON,
1385: /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON,
1386: /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM,
1387: /* GMSH_PYR */ DM_POLYTOPE_PYRAMID, DM_POLYTOPE_UNKNOWN};
1389: static inline DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType)
1390: {
1391: return DMPolytopeMap[GmshCellMap[cellType].polytope];
1392: }
1394: static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem)
1395: {
1396: DM K;
1397: PetscSpace P;
1398: PetscDualSpace Q;
1399: PetscQuadrature q, fq;
1400: PetscBool isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE;
1401: PetscBool endpoint = PETSC_TRUE;
1402: char name[32];
1404: PetscFunctionBegin;
1405: /* Create space */
1406: PetscCall(PetscSpaceCreate(comm, &P));
1407: PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL));
1408: PetscCall(PetscSpacePolynomialSetTensor(P, isTensor));
1409: PetscCall(PetscSpaceSetNumComponents(P, Nc));
1410: PetscCall(PetscSpaceSetNumVariables(P, dim));
1411: PetscCall(PetscSpaceSetDegree(P, k, PETSC_DETERMINE));
1412: if (prefix) {
1413: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix));
1414: PetscCall(PetscSpaceSetFromOptions(P));
1415: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, NULL));
1416: PetscCall(PetscSpaceGetDegree(P, &k, NULL));
1417: }
1418: PetscCall(PetscSpaceSetUp(P));
1419: /* Create dual space */
1420: PetscCall(PetscDualSpaceCreate(comm, &Q));
1421: PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE));
1422: PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor));
1423: PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity));
1424: PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0));
1425: PetscCall(PetscDualSpaceSetNumComponents(Q, Nc));
1426: PetscCall(PetscDualSpaceSetOrder(Q, k));
1427: PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K));
1428: PetscCall(PetscDualSpaceSetDM(Q, K));
1429: PetscCall(DMDestroy(&K));
1430: if (prefix) {
1431: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix));
1432: PetscCall(PetscDualSpaceSetFromOptions(Q));
1433: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, NULL));
1434: }
1435: PetscCall(PetscDualSpaceSetUp(Q));
1436: /* Create quadrature */
1437: if (isSimplex) {
1438: PetscCall(PetscDTStroudConicalQuadrature(dim, 1, k + 1, -1, +1, &q));
1439: PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, k + 1, -1, +1, &fq));
1440: } else {
1441: PetscCall(PetscDTGaussTensorQuadrature(dim, 1, k + 1, -1, +1, &q));
1442: PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, k + 1, -1, +1, &fq));
1443: }
1444: /* Create finite element */
1445: PetscCall(PetscFECreate(comm, fem));
1446: PetscCall(PetscFESetType(*fem, PETSCFEBASIC));
1447: PetscCall(PetscFESetNumComponents(*fem, Nc));
1448: PetscCall(PetscFESetBasisSpace(*fem, P));
1449: PetscCall(PetscFESetDualSpace(*fem, Q));
1450: PetscCall(PetscFESetQuadrature(*fem, q));
1451: PetscCall(PetscFESetFaceQuadrature(*fem, fq));
1452: if (prefix) {
1453: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix));
1454: PetscCall(PetscFESetFromOptions(*fem));
1455: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, NULL));
1456: }
1457: PetscCall(PetscFESetUp(*fem));
1458: /* Cleanup */
1459: PetscCall(PetscSpaceDestroy(&P));
1460: PetscCall(PetscDualSpaceDestroy(&Q));
1461: PetscCall(PetscQuadratureDestroy(&q));
1462: PetscCall(PetscQuadratureDestroy(&fq));
1463: /* Set finite element name */
1464: PetscCall(PetscSNPrintf(name, sizeof(name), "%s%" PetscInt_FMT, isSimplex ? "P" : "Q", k));
1465: PetscCall(PetscFESetName(*fem, name));
1466: PetscFunctionReturn(PETSC_SUCCESS);
1467: }
1469: /*@
1470: DMPlexCreateGmshFromFile - Create a `DMPLEX` mesh from a Gmsh file
1472: Input Parameters:
1473: + comm - The MPI communicator
1474: . filename - Name of the Gmsh file
1475: - interpolate - Create faces and edges in the mesh
1477: Output Parameter:
1478: . dm - The `DM` object representing the mesh
1480: Level: beginner
1482: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()`
1483: @*/
1484: PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1485: {
1486: PetscViewer viewer;
1487: PetscMPIInt rank;
1488: int fileType;
1489: PetscViewerType vtype;
1491: PetscFunctionBegin;
1492: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1494: /* Determine Gmsh file type (ASCII or binary) from file header */
1495: if (rank == 0) {
1496: GmshFile gmsh[1];
1497: char line[PETSC_MAX_PATH_LEN];
1498: int snum;
1499: float version;
1500: int fileFormat;
1502: PetscCall(PetscArrayzero(gmsh, 1));
1503: PetscCall(PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer));
1504: PetscCall(PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII));
1505: PetscCall(PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ));
1506: PetscCall(PetscViewerFileSetName(gmsh->viewer, filename));
1507: /* Read only the first two lines of the Gmsh file */
1508: PetscCall(GmshReadSection(gmsh, line));
1509: PetscCall(GmshExpect(gmsh, "$MeshFormat", line));
1510: PetscCall(GmshReadString(gmsh, line, 2));
1511: snum = sscanf(line, "%f %d", &version, &fileType);
1512: fileFormat = (int)roundf(version * 10);
1513: PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line);
1514: PetscCheck(fileFormat >= 22, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version);
1515: PetscCheck((int)version != 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version);
1516: PetscCheck(fileFormat <= 41, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version);
1517: PetscCall(PetscViewerDestroy(&gmsh->viewer));
1518: }
1519: PetscCallMPI(MPI_Bcast(&fileType, 1, MPI_INT, 0, comm));
1520: vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY;
1522: /* Create appropriate viewer and build plex */
1523: PetscCall(PetscViewerCreate(comm, &viewer));
1524: PetscCall(PetscViewerSetType(viewer, vtype));
1525: PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
1526: PetscCall(PetscViewerFileSetName(viewer, filename));
1527: PetscCall(DMPlexCreateGmsh(comm, viewer, interpolate, dm));
1528: PetscCall(PetscViewerDestroy(&viewer));
1529: PetscFunctionReturn(PETSC_SUCCESS);
1530: }
1532: /*@
1533: DMPlexCreateGmsh - Create a `DMPLEX` mesh from a Gmsh file viewer
1535: Collective
1537: Input Parameters:
1538: + comm - The MPI communicator
1539: . viewer - The `PetscViewer` associated with a Gmsh file
1540: - interpolate - Create faces and edges in the mesh
1542: Output Parameter:
1543: . dm - The `DM` object representing the mesh
1545: Options Database Keys:
1546: + -dm_plex_gmsh_hybrid - Force triangular prisms to use tensor order
1547: . -dm_plex_gmsh_periodic - Read Gmsh periodic section and construct a periodic Plex
1548: . -dm_plex_gmsh_highorder - Generate high-order coordinates
1549: . -dm_plex_gmsh_project - Project high-order coordinates to a different space, use the prefix dm_plex_gmsh_project_ to define the space
1550: . -dm_plex_gmsh_use_generic - Generate generic labels, i.e. Cell Sets, Face Sets, etc.
1551: . -dm_plex_gmsh_use_regions - Generate labels with region names
1552: . -dm_plex_gmsh_mark_vertices - Add vertices to generated labels
1553: . -dm_plex_gmsh_mark_vertices_strict - Add vertices included in a region to generated labels
1554: . -dm_plex_gmsh_multiple_tags - Allow multiple tags for default labels
1555: - -dm_plex_gmsh_spacedim d - Embedding space dimension, if different from topological dimension
1557: Level: beginner
1559: Notes:
1560: The Gmsh file format is described in <http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format>
1562: By default, the "Cell Sets", "Face Sets", and "Vertex Sets" labels are created, and only insert the first tag on a point. By using `-dm_plex_gmsh_multiple_tags`, all tags can be inserted. Alternatively, `-dm_plex_gmsh_use_regions` creates labels based on the region names from the `PhysicalNames` section, and all tags are used, but you can retain the generic labels using `-dm_plex_gmsh_use_generic`.
1564: .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`
1565: @*/
1566: PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
1567: {
1568: GmshMesh *mesh = NULL;
1569: PetscViewer parentviewer = NULL;
1570: PetscBT periodicVerts = NULL;
1571: PetscBT *periodicCells = NULL;
1572: DM cdm, cdmCell = NULL;
1573: PetscSection cs, csCell = NULL;
1574: Vec coordinates, coordinatesCell;
1575: DMLabel cellSets = NULL, faceSets = NULL, edgeSets = NULL, vertSets = NULL, marker = NULL, *regionSets;
1576: PetscInt dim = 0, coordDim = -1, order = 0, maxHeight = 0;
1577: PetscInt numNodes = 0, numElems = 0, numVerts = 0, numCells = 0, vStart, vEnd;
1578: PetscInt cell, cone[8], e, n, v, d;
1579: PetscBool usegeneric = PETSC_TRUE, useregions = PETSC_FALSE, markvertices = PETSC_FALSE, markverticesstrict = PETSC_FALSE, multipleTags = PETSC_FALSE;
1580: PetscBool flg, binary, hybrid = interpolate, periodic = PETSC_TRUE;
1581: PetscBool highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE;
1582: PetscBool isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE;
1583: PetscMPIInt rank;
1585: PetscFunctionBegin;
1586: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1587: PetscObjectOptionsBegin((PetscObject)viewer);
1588: PetscOptionsHeadBegin(PetscOptionsObject, "DMPlex Gmsh options");
1589: PetscCall(PetscOptionsBool("-dm_plex_gmsh_hybrid", "Force triangular prisms to use tensor order", "DMPlexCreateGmsh", hybrid, &hybrid, NULL));
1590: PetscCall(PetscOptionsBool("-dm_plex_gmsh_periodic", "Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL));
1591: PetscCall(PetscOptionsBool("-dm_plex_gmsh_highorder", "Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet));
1592: PetscCall(PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL));
1593: PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_regions", "Generate labels with region names", "DMPlexCreateGmsh", useregions, &useregions, NULL));
1594: PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_generic", "Generate generic labels, i.e. Cell Sets, Face Sets, etc", "DMPlexCreateGmsh", usegeneric, &usegeneric, &flg));
1595: if (!flg && useregions) usegeneric = PETSC_FALSE;
1596: PetscCall(PetscOptionsBool("-dm_plex_gmsh_mark_vertices", "Add vertices to generated labels", "DMPlexCreateGmsh", markvertices, &markvertices, NULL));
1597: PetscCall(PetscOptionsBool("-dm_plex_gmsh_mark_vertices_strict", "Add only directly tagged vertices to generated labels", "DMPlexCreateGmsh", markverticesstrict, &markverticesstrict, NULL));
1598: PetscCall(PetscOptionsBool("-dm_plex_gmsh_multiple_tags", "Allow multiple tags for default labels", "DMPlexCreateGmsh", multipleTags, &multipleTags, NULL));
1599: PetscCall(PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE));
1600: PetscCall(PetscOptionsBoundedInt("-dm_localize_height", "Localize edges and faces in addition to cells", "", maxHeight, &maxHeight, NULL, 0));
1601: PetscOptionsHeadEnd();
1602: PetscOptionsEnd();
1604: PetscCall(GmshCellInfoSetUp());
1606: PetscCall(DMCreate(comm, dm));
1607: PetscCall(DMSetType(*dm, DMPLEX));
1608: PetscCall(PetscLogEventBegin(DMPLEX_CreateGmsh, *dm, NULL, NULL, NULL));
1610: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary));
1612: /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */
1613: if (binary) {
1614: parentviewer = viewer;
1615: PetscCall(PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer));
1616: }
1618: if (rank == 0) {
1619: GmshFile gmsh[1];
1620: char line[PETSC_MAX_PATH_LEN];
1621: PetscBool match;
1623: PetscCall(PetscArrayzero(gmsh, 1));
1624: gmsh->viewer = viewer;
1625: gmsh->binary = binary;
1627: PetscCall(GmshMeshCreate(&mesh));
1629: /* Read mesh format */
1630: PetscCall(GmshReadSection(gmsh, line));
1631: PetscCall(GmshExpect(gmsh, "$MeshFormat", line));
1632: PetscCall(GmshReadMeshFormat(gmsh));
1633: PetscCall(GmshReadEndSection(gmsh, "$EndMeshFormat", line));
1635: /* OPTIONAL Read mesh version (Neper only) */
1636: PetscCall(GmshReadSection(gmsh, line));
1637: PetscCall(GmshMatch(gmsh, "$MeshVersion", line, &match));
1638: if (match) {
1639: PetscCall(GmshExpect(gmsh, "$MeshVersion", line));
1640: PetscCall(GmshReadMeshVersion(gmsh));
1641: PetscCall(GmshReadEndSection(gmsh, "$EndMeshVersion", line));
1642: /* Initial read for entity section */
1643: PetscCall(GmshReadSection(gmsh, line));
1644: }
1646: /* OPTIONAL Read mesh domain (Neper only) */
1647: PetscCall(GmshMatch(gmsh, "$Domain", line, &match));
1648: if (match) {
1649: PetscCall(GmshExpect(gmsh, "$Domain", line));
1650: PetscCall(GmshReadMeshDomain(gmsh));
1651: PetscCall(GmshReadEndSection(gmsh, "$EndDomain", line));
1652: /* Initial read for entity section */
1653: PetscCall(GmshReadSection(gmsh, line));
1654: }
1656: /* OPTIONAL Read physical names */
1657: PetscCall(GmshMatch(gmsh, "$PhysicalNames", line, &match));
1658: if (match) {
1659: PetscCall(GmshExpect(gmsh, "$PhysicalNames", line));
1660: PetscCall(GmshReadPhysicalNames(gmsh, mesh));
1661: PetscCall(GmshReadEndSection(gmsh, "$EndPhysicalNames", line));
1662: /* Initial read for entity section */
1663: PetscCall(GmshReadSection(gmsh, line));
1664: }
1666: /* OPTIONAL Read entities */
1667: if (gmsh->fileFormat >= 40) {
1668: PetscCall(GmshMatch(gmsh, "$Entities", line, &match));
1669: if (match) {
1670: PetscCall(GmshExpect(gmsh, "$Entities", line));
1671: PetscCall(GmshReadEntities(gmsh, mesh));
1672: PetscCall(GmshReadEndSection(gmsh, "$EndEntities", line));
1673: /* Initial read for nodes section */
1674: PetscCall(GmshReadSection(gmsh, line));
1675: }
1676: }
1678: /* Read nodes */
1679: PetscCall(GmshExpect(gmsh, "$Nodes", line));
1680: PetscCall(GmshReadNodes(gmsh, mesh));
1681: PetscCall(GmshReadEndSection(gmsh, "$EndNodes", line));
1683: /* Read elements */
1684: PetscCall(GmshReadSection(gmsh, line));
1685: PetscCall(GmshExpect(gmsh, "$Elements", line));
1686: PetscCall(GmshReadElements(gmsh, mesh));
1687: PetscCall(GmshReadEndSection(gmsh, "$EndElements", line));
1689: /* Read periodic section (OPTIONAL) */
1690: if (periodic) {
1691: PetscCall(GmshReadSection(gmsh, line));
1692: PetscCall(GmshMatch(gmsh, "$Periodic", line, &periodic));
1693: }
1694: if (periodic) {
1695: PetscCall(GmshExpect(gmsh, "$Periodic", line));
1696: PetscCall(GmshReadPeriodic(gmsh, mesh));
1697: PetscCall(GmshReadEndSection(gmsh, "$EndPeriodic", line));
1698: }
1700: PetscCall(PetscFree(gmsh->wbuf));
1701: PetscCall(PetscFree(gmsh->sbuf));
1702: PetscCall(PetscFree(gmsh->nbuf));
1704: dim = mesh->dim;
1705: order = mesh->order;
1706: numNodes = mesh->numNodes;
1707: numElems = mesh->numElems;
1708: numVerts = mesh->numVerts;
1709: numCells = mesh->numCells;
1711: {
1712: GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL;
1713: GmshElement *elemB = PetscSafePointerPlusOffset(elemA, mesh->numCells - 1);
1714: int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1;
1715: int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1;
1716: isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE;
1717: isHybrid = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE;
1718: hasTetra = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE;
1719: }
1720: }
1722: if (parentviewer) PetscCall(PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer));
1724: {
1725: int buf[6];
1727: buf[0] = (int)dim;
1728: buf[1] = (int)order;
1729: buf[2] = periodic;
1730: buf[3] = isSimplex;
1731: buf[4] = isHybrid;
1732: buf[5] = hasTetra;
1734: PetscCallMPI(MPI_Bcast(buf, 6, MPI_INT, 0, comm));
1736: dim = buf[0];
1737: order = buf[1];
1738: periodic = buf[2] ? PETSC_TRUE : PETSC_FALSE;
1739: isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE;
1740: isHybrid = buf[4] ? PETSC_TRUE : PETSC_FALSE;
1741: hasTetra = buf[5] ? PETSC_TRUE : PETSC_FALSE;
1742: }
1744: if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE;
1745: PetscCheck(!highOrder || !isHybrid, comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet");
1747: /* We do not want this label automatically computed, instead we fill it here */
1748: PetscCall(DMCreateLabel(*dm, "celltype"));
1750: /* Allocate the cell-vertex mesh */
1751: PetscCall(DMPlexSetChart(*dm, 0, numCells + numVerts));
1752: for (cell = 0; cell < numCells; ++cell) {
1753: GmshElement *elem = mesh->elements + cell;
1754: DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType);
1755: if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR;
1756: PetscCall(DMPlexSetConeSize(*dm, cell, elem->numVerts));
1757: PetscCall(DMPlexSetCellType(*dm, cell, ctype));
1758: }
1759: for (v = numCells; v < numCells + numVerts; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
1760: PetscCall(DMSetUp(*dm));
1762: /* Add cell-vertex connections */
1763: for (cell = 0; cell < numCells; ++cell) {
1764: GmshElement *elem = mesh->elements + cell;
1765: for (v = 0; v < elem->numVerts; ++v) {
1766: const PetscInt nn = elem->nodes[v];
1767: const PetscInt vv = mesh->vertexMap[nn];
1768: cone[v] = numCells + vv;
1769: }
1770: PetscCall(DMPlexReorderCell(*dm, cell, cone));
1771: PetscCall(DMPlexSetCone(*dm, cell, cone));
1772: }
1774: PetscCall(DMSetDimension(*dm, dim));
1775: PetscCall(DMPlexSymmetrize(*dm));
1776: PetscCall(DMPlexStratify(*dm));
1777: if (interpolate) {
1778: DM idm;
1780: PetscCall(DMPlexInterpolate(*dm, &idm));
1781: PetscCall(DMDestroy(dm));
1782: *dm = idm;
1783: }
1784: PetscCall(DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd));
1786: if (rank == 0) {
1787: const PetscInt Nr = useregions ? mesh->numRegions : 0;
1789: PetscCall(PetscCalloc1(Nr, ®ionSets));
1790: for (cell = 0, e = 0; e < numElems; ++e) {
1791: GmshElement *elem = mesh->elements + e;
1793: /* Create cell sets */
1794: if (elem->dim == dim && dim > 0) {
1795: if (elem->numTags > 0) {
1796: PetscInt Nt = elem->numTags, t, r;
1798: for (t = 0; t < Nt; ++t) {
1799: const PetscInt tag = elem->tags[t];
1800: const PetscBool generic = usegeneric && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE;
1802: if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", cell, tag));
1803: for (r = 0; r < Nr; ++r) {
1804: if (mesh->regionDims[r] != dim) continue;
1805: if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], cell, tag));
1806: }
1807: }
1808: }
1809: cell++;
1810: }
1812: /* Create face sets */
1813: if (elem->numTags && interpolate && elem->dim == dim - 1) {
1814: PetscInt joinSize;
1815: const PetscInt *join = NULL;
1816: PetscInt Nt = elem->numTags, pdepth;
1818: /* Find the relevant facet with vertex joins */
1819: for (v = 0; v < elem->numVerts; ++v) {
1820: const PetscInt nn = elem->nodes[v];
1821: const PetscInt vv = mesh->vertexMap[nn];
1822: cone[v] = vStart + vv;
1823: }
1824: PetscCall(DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join));
1825: PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Could not determine Plex facet for Gmsh element %" PetscInt_FMT " (Plex cell %" PetscInt_FMT ")", elem->id, e);
1826: PetscCall(DMPlexGetPointDepth(*dm, join[0], &pdepth));
1827: PetscCheck(pdepth == dim - 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Plex facet %" PetscInt_FMT " for Gmsh element %" PetscInt_FMT " had depth %" PetscInt_FMT " != %" PetscInt_FMT, join[0], elem->id, pdepth, dim - 1);
1828: for (PetscInt t = 0; t < Nt; ++t) {
1829: const PetscInt tag = elem->tags[t];
1830: const PetscBool generic = usegeneric && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE;
1832: if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], tag));
1833: for (PetscInt r = 0; r < Nr; ++r) {
1834: if (mesh->regionDims[r] != dim - 1) continue;
1835: if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], join[0], tag));
1836: }
1837: }
1838: PetscCall(DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join));
1839: }
1841: /* Create edge sets */
1842: if (elem->numTags && interpolate && dim > 2 && elem->dim == 1) {
1843: PetscInt joinSize;
1844: const PetscInt *join = NULL;
1845: PetscInt Nt = elem->numTags, t, r;
1847: /* Find the relevant edge with vertex joins */
1848: for (v = 0; v < elem->numVerts; ++v) {
1849: const PetscInt nn = elem->nodes[v];
1850: const PetscInt vv = mesh->vertexMap[nn];
1851: cone[v] = vStart + vv;
1852: }
1853: PetscCall(DMPlexGetJoin(*dm, elem->numVerts, cone, &joinSize, &join));
1854: PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Could not determine Plex edge for Gmsh element %" PetscInt_FMT " (Plex cell %" PetscInt_FMT ")", elem->id, e);
1855: for (t = 0; t < Nt; ++t) {
1856: const PetscInt tag = elem->tags[t];
1857: const PetscBool generic = usegeneric && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE;
1859: if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &edgeSets, "Edge Sets", join[0], tag));
1860: for (r = 0; r < Nr; ++r) {
1861: if (mesh->regionDims[r] != 1) continue;
1862: if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], join[0], tag));
1863: }
1864: }
1865: PetscCall(DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join));
1866: }
1868: /* Create vertex sets */
1869: if (elem->numTags && elem->dim == 0 && (markverticesstrict || markvertices)) {
1870: const PetscInt nn = elem->nodes[0];
1871: const PetscInt vv = mesh->vertexMap[nn];
1872: PetscInt Nt = elem->numTags;
1874: for (PetscInt t = 0; t < Nt; ++t) {
1875: const PetscInt tag = elem->tags[t];
1877: if (vv < 0) continue;
1878: if (usegeneric) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag));
1879: for (PetscInt r = 0; r < Nr; ++r) {
1880: if (mesh->regionDims[r] != 0) continue;
1881: if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag));
1882: }
1883: }
1884: }
1885: }
1886: if (markvertices) {
1887: for (v = 0; v < numNodes; ++v) {
1888: const PetscInt vv = mesh->vertexMap[v];
1889: const PetscInt *tags = &mesh->nodelist->tag[v * GMSH_MAX_TAGS];
1891: if (vv < 0) continue;
1892: for (PetscInt t = 0; t < GMSH_MAX_TAGS; ++t) {
1893: const PetscInt tag = tags[t];
1894: const PetscBool generic = usegeneric && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE;
1896: if (tag == -1) continue;
1897: if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag));
1898: for (PetscInt r = 0; r < Nr; ++r) {
1899: if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag));
1900: }
1901: }
1902: }
1903: }
1904: PetscCall(PetscFree(regionSets));
1905: }
1907: { /* Create Cell/Face/Edge/Vertex Sets labels at all processes */
1908: enum {
1909: n = 5
1910: };
1911: PetscBool flag[n];
1913: flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1914: flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1915: flag[2] = edgeSets ? PETSC_TRUE : PETSC_FALSE;
1916: flag[3] = vertSets ? PETSC_TRUE : PETSC_FALSE;
1917: flag[4] = marker ? PETSC_TRUE : PETSC_FALSE;
1918: PetscCallMPI(MPI_Bcast(flag, n, MPI_C_BOOL, 0, comm));
1919: if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
1920: if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
1921: if (flag[2]) PetscCall(DMCreateLabel(*dm, "Edge Sets"));
1922: if (flag[3]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
1923: if (flag[4]) PetscCall(DMCreateLabel(*dm, "marker"));
1924: }
1926: if (periodic) {
1927: PetscCall(PetscBTCreate(numVerts, &periodicVerts));
1928: for (n = 0; n < numNodes; ++n) {
1929: if (mesh->vertexMap[n] >= 0) {
1930: if (PetscUnlikely(mesh->periodMap[n] != n)) {
1931: PetscInt m = mesh->periodMap[n];
1932: PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[n]));
1933: PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[m]));
1934: }
1935: }
1936: }
1937: PetscCall(DMGetCoordinateDM(*dm, &cdm));
1938: PetscCall(PetscMalloc1(maxHeight + 1, &periodicCells));
1939: for (PetscInt h = 0; h <= maxHeight; ++h) {
1940: PetscInt pStart, pEnd;
1942: PetscCall(DMPlexGetHeightStratum(*dm, h, &pStart, &pEnd));
1943: PetscCall(PetscBTCreate(pEnd - pStart, &periodicCells[h]));
1944: for (PetscInt p = pStart; p < pEnd; ++p) {
1945: PetscInt *closure = NULL;
1946: PetscInt Ncl;
1948: PetscCall(DMPlexGetTransitiveClosure(*dm, p, PETSC_TRUE, &Ncl, &closure));
1949: for (PetscInt cl = 0; cl < Ncl * 2; cl += 2) {
1950: if (closure[cl] >= vStart && closure[cl] < vEnd) {
1951: if (PetscUnlikely(PetscBTLookup(periodicVerts, closure[cl] - vStart))) {
1952: PetscCall(PetscBTSet(periodicCells[h], p - pStart));
1953: break;
1954: }
1955: }
1956: }
1957: PetscCall(DMPlexRestoreTransitiveClosure(*dm, p, PETSC_TRUE, &Ncl, &closure));
1958: }
1959: }
1960: }
1962: /* Setup coordinate DM */
1963: if (coordDim < 0) coordDim = dim;
1964: PetscCall(DMSetCoordinateDim(*dm, coordDim));
1965: PetscCall(DMGetCoordinateDM(*dm, &cdm));
1966: if (highOrder) {
1967: PetscFE fe;
1968: PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
1969: PetscDTNodeType nodeType = PETSCDTNODES_EQUISPACED;
1971: PetscCall(GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe));
1972: PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view"));
1973: PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)fe));
1974: PetscCall(PetscFEDestroy(&fe));
1975: PetscCall(DMCreateDS(cdm));
1976: }
1977: if (periodic) {
1978: PetscCall(DMClone(cdm, &cdmCell));
1979: PetscCall(DMSetCellCoordinateDM(*dm, cdmCell));
1980: }
1982: /* Create coordinates */
1983: if (highOrder) {
1984: PetscInt maxDof = GmshNumNodes_HEX(order) * coordDim;
1985: double *coords = mesh ? mesh->nodelist->xyz : NULL;
1986: PetscSection section;
1987: PetscScalar *cellCoords;
1989: PetscCall(DMSetLocalSection(cdm, NULL));
1990: PetscCall(DMGetLocalSection(cdm, &cs));
1991: PetscCall(PetscSectionClone(cs, §ion));
1992: if (isSimplex) PetscCall(DMPlexSetClosurePermutationLexicographic(cdm, 0, section));
1993: else PetscCall(DMPlexSetClosurePermutationTensor(cdm, 0, section));
1995: PetscCall(DMCreateLocalVector(cdm, &coordinates));
1996: PetscCall(PetscMalloc1(maxDof, &cellCoords));
1997: for (cell = 0; cell < numCells; ++cell) {
1998: GmshElement *elem = mesh->elements + cell;
1999: const int *lexorder = GmshCellMap[elem->cellType].lexorder();
2000: int s = 0;
2001: for (n = 0; n < elem->numNodes; ++n) {
2002: while (lexorder[n + s] < 0) ++s;
2003: const PetscInt node = elem->nodes[lexorder[n + s]];
2004: for (d = 0; d < coordDim; ++d) cellCoords[(n + s) * coordDim + d] = (PetscReal)coords[node * 3 + d];
2005: }
2006: if (s) {
2007: /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */
2008: PetscReal quaCenterWeights[9] = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25};
2009: /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */
2010: PetscReal hexBottomWeights[27] = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
2011: PetscReal hexFrontWeights[27] = {-0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0};
2012: PetscReal hexLeftWeights[27] = {-0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0};
2013: PetscReal hexRightWeights[27] = {0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25};
2014: PetscReal hexBackWeights[27] = {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25};
2015: PetscReal hexTopWeights[27] = {0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25};
2016: PetscReal hexCenterWeights[27] = {-0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25, 0.25, 0.0, 0.25, 0.0, 0.0, 0.0, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25};
2017: PetscReal *sdWeights2[9] = {NULL, NULL, NULL, NULL, quaCenterWeights, NULL, NULL, NULL, NULL};
2018: PetscReal *sdWeights3[27] = {NULL, NULL, NULL, NULL, hexBottomWeights, NULL, NULL, NULL, NULL, NULL, hexFrontWeights, NULL, hexLeftWeights, hexCenterWeights, hexRightWeights, NULL, hexBackWeights, NULL,
2019: NULL, NULL, NULL, NULL, hexTopWeights, NULL, NULL, NULL, NULL};
2020: PetscReal **sdWeights[4] = {NULL, NULL, sdWeights2, sdWeights3};
2022: /* Missing entries in serendipity cell, only works for 8-node quad and 20-node hex */
2023: for (n = 0; n < elem->numNodes + s; ++n) {
2024: if (lexorder[n] >= 0) continue;
2025: for (d = 0; d < coordDim; ++d) cellCoords[n * coordDim + d] = 0.0;
2026: for (int bn = 0; bn < elem->numNodes + s; ++bn) {
2027: if (lexorder[bn] < 0) continue;
2028: const PetscReal *weights = sdWeights[coordDim][n];
2029: const PetscInt bnode = elem->nodes[lexorder[bn]];
2030: for (d = 0; d < coordDim; ++d) cellCoords[n * coordDim + d] += weights[bn] * (PetscReal)coords[bnode * 3 + d];
2031: }
2032: }
2033: }
2034: PetscCall(DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES));
2035: }
2036: PetscCall(PetscSectionDestroy(§ion));
2037: PetscCall(PetscFree(cellCoords));
2038: } else {
2039: PetscInt *nodeMap;
2040: double *coords = mesh ? mesh->nodelist->xyz : NULL;
2041: PetscScalar *pointCoords;
2043: PetscCall(DMGetCoordinateSection(*dm, &cs));
2044: PetscCall(PetscSectionSetNumFields(cs, 1));
2045: PetscCall(PetscSectionSetFieldComponents(cs, 0, coordDim));
2046: PetscCall(PetscSectionSetChart(cs, numCells, numCells + numVerts));
2047: for (v = numCells; v < numCells + numVerts; ++v) {
2048: PetscCall(PetscSectionSetDof(cs, v, coordDim));
2049: PetscCall(PetscSectionSetFieldDof(cs, v, 0, coordDim));
2050: }
2051: PetscCall(PetscSectionSetUp(cs));
2053: /* We need to localize coordinates on cells */
2054: if (periodic) {
2055: PetscInt newStart = PETSC_INT_MAX, newEnd = -1, pStart, pEnd;
2057: PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)cdmCell), &csCell));
2058: PetscCall(PetscSectionSetNumFields(csCell, 1));
2059: PetscCall(PetscSectionSetFieldComponents(csCell, 0, coordDim));
2060: for (PetscInt h = 0; h <= maxHeight; h++) {
2061: PetscCall(DMPlexGetHeightStratum(cdmCell, h, &pStart, &pEnd));
2062: newStart = PetscMin(newStart, pStart);
2063: newEnd = PetscMax(newEnd, pEnd);
2064: }
2065: PetscCall(PetscSectionSetChart(csCell, newStart, newEnd));
2066: for (PetscInt h = 0; h <= maxHeight; h++) {
2067: PetscCall(DMPlexGetHeightStratum(cdmCell, h, &pStart, &pEnd));
2068: for (PetscInt p = pStart; p < pEnd; ++p) {
2069: PetscInt *closure = NULL;
2070: PetscInt Ncl, Nv = 0;
2072: if (PetscUnlikely(PetscBTLookup(periodicCells[h], p - pStart))) {
2073: PetscCall(DMPlexGetTransitiveClosure(*dm, p, PETSC_TRUE, &Ncl, &closure));
2074: for (PetscInt cl = 0; cl < Ncl * 2; cl += 2) {
2075: if (closure[cl] >= vStart && closure[cl] < vEnd) ++Nv;
2076: }
2077: PetscCall(DMPlexRestoreTransitiveClosure(*dm, p, PETSC_TRUE, &Ncl, &closure));
2078: PetscCall(PetscSectionSetDof(csCell, p, Nv * coordDim));
2079: PetscCall(PetscSectionSetFieldDof(csCell, p, 0, Nv * coordDim));
2080: }
2081: }
2082: }
2083: PetscCall(PetscSectionSetUp(csCell));
2084: PetscCall(DMSetCellCoordinateSection(*dm, PETSC_DETERMINE, csCell));
2085: }
2087: PetscCall(DMCreateLocalVector(cdm, &coordinates));
2088: PetscCall(VecGetArray(coordinates, &pointCoords));
2089: PetscCall(PetscMalloc1(numVerts, &nodeMap));
2090: for (n = 0; n < numNodes; n++)
2091: if (mesh->vertexMap[n] >= 0) nodeMap[mesh->vertexMap[n]] = n;
2092: for (v = 0; v < numVerts; ++v) {
2093: PetscInt off, node = nodeMap[v];
2095: PetscCall(PetscSectionGetOffset(cs, numCells + v, &off));
2096: for (d = 0; d < coordDim; ++d) pointCoords[off + d] = (PetscReal)coords[node * 3 + d];
2097: }
2098: PetscCall(VecRestoreArray(coordinates, &pointCoords));
2100: if (periodic) {
2101: PetscInt cStart, cEnd;
2103: PetscCall(DMPlexGetHeightStratum(cdmCell, 0, &cStart, &cEnd));
2104: PetscCall(DMCreateLocalVector(cdmCell, &coordinatesCell));
2105: PetscCall(VecGetArray(coordinatesCell, &pointCoords));
2106: for (PetscInt c = cStart; c < cEnd; ++c) {
2107: GmshElement *elem = mesh->elements + c - cStart;
2108: PetscInt *closure = NULL;
2109: PetscInt verts[8];
2110: PetscInt Ncl, Nv = 0;
2112: for (PetscInt v = 0; v < elem->numVerts; ++v) cone[v] = elem->nodes[v];
2113: PetscCall(DMPlexReorderCell(cdmCell, c, cone));
2114: PetscCall(DMPlexGetTransitiveClosure(cdmCell, c, PETSC_TRUE, &Ncl, &closure));
2115: for (PetscInt cl = 0; cl < Ncl * 2; cl += 2) {
2116: if (closure[cl] >= vStart && closure[cl] < vEnd) verts[Nv++] = closure[cl];
2117: }
2118: PetscCheck(Nv == elem->numVerts, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of vertices %" PetscInt_FMT " in closure does not match number of vertices %" PetscInt_FMT " in Gmsh cell", Nv, elem->numVerts);
2119: for (PetscInt cl = 0; cl < Ncl * 2; cl += 2) {
2120: const PetscInt point = closure[cl];
2121: PetscInt hStart, h;
2123: PetscCall(DMPlexGetPointHeight(cdmCell, point, &h));
2124: if (h > maxHeight) continue;
2125: PetscCall(DMPlexGetHeightStratum(cdmCell, h, &hStart, NULL));
2126: if (PetscUnlikely(PetscBTLookup(periodicCells[h], point - hStart))) {
2127: PetscInt *pclosure = NULL;
2128: PetscInt Npcl, off, v;
2130: PetscCall(PetscSectionGetOffset(csCell, point, &off));
2131: PetscCall(DMPlexGetTransitiveClosure(cdmCell, point, PETSC_TRUE, &Npcl, &pclosure));
2132: for (PetscInt pcl = 0; pcl < Npcl * 2; pcl += 2) {
2133: if (pclosure[pcl] >= vStart && pclosure[pcl] < vEnd) {
2134: for (v = 0; v < Nv; ++v)
2135: if (verts[v] == pclosure[pcl]) break;
2136: PetscCheck(v < Nv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not find vertex %" PetscInt_FMT " in closure of cell %" PetscInt_FMT, pclosure[pcl], c);
2137: for (PetscInt d = 0; d < coordDim; ++d) pointCoords[off++] = (PetscReal)coords[cone[v] * 3 + d];
2138: ++v;
2139: }
2140: }
2141: PetscCall(DMPlexRestoreTransitiveClosure(cdmCell, point, PETSC_TRUE, &Npcl, &pclosure));
2142: }
2143: }
2144: PetscCall(DMPlexRestoreTransitiveClosure(cdmCell, c, PETSC_TRUE, &Ncl, &closure));
2145: }
2146: PetscCall(VecSetBlockSize(coordinatesCell, coordDim));
2147: PetscCall(VecRestoreArray(coordinatesCell, &pointCoords));
2148: PetscCall(DMSetCellCoordinatesLocal(*dm, coordinatesCell));
2149: PetscCall(VecDestroy(&coordinatesCell));
2150: }
2151: PetscCall(PetscFree(nodeMap));
2152: PetscCall(PetscSectionDestroy(&csCell));
2153: PetscCall(DMDestroy(&cdmCell));
2154: }
2156: PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
2157: PetscCall(VecSetBlockSize(coordinates, coordDim));
2158: PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
2159: PetscCall(VecDestroy(&coordinates));
2161: PetscCall(GmshMeshDestroy(&mesh));
2162: PetscCall(PetscBTDestroy(&periodicVerts));
2163: if (periodic) {
2164: for (PetscInt h = 0; h <= maxHeight; ++h) PetscCall(PetscBTDestroy(periodicCells + h));
2165: PetscCall(PetscFree(periodicCells));
2166: }
2168: if (highOrder && project) {
2169: PetscFE fe;
2170: const char prefix[] = "dm_plex_gmsh_project_";
2171: PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE;
2172: PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI;
2174: PetscCall(GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe));
2175: PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view"));
2176: PetscCall(DMSetCoordinateDisc(*dm, fe, PETSC_FALSE, PETSC_TRUE));
2177: PetscCall(PetscFEDestroy(&fe));
2178: }
2180: PetscCall(PetscLogEventEnd(DMPLEX_CreateGmsh, *dm, NULL, NULL, NULL));
2181: PetscFunctionReturn(PETSC_SUCCESS);
2182: }