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) \ 89371c9d4SSatish Balay static int *GmshLexOrder_##T##_##p(void) { \ 9066ea43fSLisandro Dalcin static int Gmsh_LexOrder_##T##_##p[GmshNumNodes_##T(p)] = {-1}; \ 10066ea43fSLisandro Dalcin int *lex = Gmsh_LexOrder_##T##_##p; \ 11066ea43fSLisandro Dalcin if (lex[0] == -1) (void)GmshLexOrder_##T(p, lex, 0); \ 12066ea43fSLisandro Dalcin return lex; \ 13066ea43fSLisandro Dalcin } 14066ea43fSLisandro Dalcin 159371c9d4SSatish Balay static int *GmshLexOrder_QUA_2_Serendipity(void) { 16b9bf55e5SMatthew G. Knepley static int Gmsh_LexOrder_QUA_2_Serendipity[9] = {-1}; 17b9bf55e5SMatthew G. Knepley int *lex = Gmsh_LexOrder_QUA_2_Serendipity; 18b9bf55e5SMatthew G. Knepley if (lex[0] == -1) { 19b9bf55e5SMatthew G. Knepley /* Vertices */ 209371c9d4SSatish Balay lex[0] = 0; 219371c9d4SSatish Balay lex[2] = 1; 229371c9d4SSatish Balay lex[8] = 2; 239371c9d4SSatish Balay lex[6] = 3; 24b9bf55e5SMatthew G. Knepley /* Edges */ 259371c9d4SSatish Balay lex[1] = 4; 269371c9d4SSatish Balay lex[5] = 5; 279371c9d4SSatish Balay lex[7] = 6; 289371c9d4SSatish Balay lex[3] = 7; 29b9bf55e5SMatthew G. Knepley /* Cell */ 30b9bf55e5SMatthew G. Knepley lex[4] = -8; 31b9bf55e5SMatthew G. Knepley } 32b9bf55e5SMatthew G. Knepley return lex; 33b9bf55e5SMatthew G. Knepley } 34b9bf55e5SMatthew G. Knepley 359371c9d4SSatish Balay static int *GmshLexOrder_HEX_2_Serendipity(void) { 36b9bf55e5SMatthew G. Knepley static int Gmsh_LexOrder_HEX_2_Serendipity[27] = {-1}; 37b9bf55e5SMatthew G. Knepley int *lex = Gmsh_LexOrder_HEX_2_Serendipity; 38b9bf55e5SMatthew G. Knepley if (lex[0] == -1) { 39b9bf55e5SMatthew G. Knepley /* Vertices */ 409371c9d4SSatish Balay lex[0] = 0; 419371c9d4SSatish Balay lex[2] = 1; 429371c9d4SSatish Balay lex[8] = 2; 439371c9d4SSatish Balay lex[6] = 3; 449371c9d4SSatish Balay lex[18] = 4; 459371c9d4SSatish Balay lex[20] = 5; 469371c9d4SSatish Balay lex[26] = 6; 479371c9d4SSatish Balay lex[24] = 7; 48b9bf55e5SMatthew G. Knepley /* Edges */ 499371c9d4SSatish Balay lex[1] = 8; 509371c9d4SSatish Balay lex[3] = 9; 519371c9d4SSatish Balay lex[9] = 10; 529371c9d4SSatish Balay lex[5] = 11; 539371c9d4SSatish Balay lex[11] = 12; 549371c9d4SSatish Balay lex[7] = 13; 559371c9d4SSatish Balay lex[17] = 14; 569371c9d4SSatish Balay lex[15] = 15; 579371c9d4SSatish Balay lex[19] = 16; 589371c9d4SSatish Balay lex[21] = 17; 599371c9d4SSatish Balay lex[23] = 18; 609371c9d4SSatish Balay lex[25] = 19; 61b9bf55e5SMatthew G. Knepley /* Faces */ 629371c9d4SSatish Balay lex[4] = -20; 639371c9d4SSatish Balay lex[10] = -21; 649371c9d4SSatish Balay lex[12] = -22; 659371c9d4SSatish Balay lex[14] = -23; 669371c9d4SSatish Balay lex[16] = -24; 679371c9d4SSatish Balay lex[22] = -25; 68b9bf55e5SMatthew G. Knepley /* Cell */ 69b9bf55e5SMatthew G. Knepley lex[13] = -26; 70b9bf55e5SMatthew G. Knepley } 71b9bf55e5SMatthew G. Knepley return lex; 72b9bf55e5SMatthew G. Knepley } 73b9bf55e5SMatthew G. Knepley 74066ea43fSLisandro Dalcin #define GMSH_LEXORDER_LIST(T) \ 75066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 1) \ 76066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 2) \ 77066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 3) \ 78066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 4) \ 79066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 5) \ 80066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 6) \ 81066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 7) \ 82066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 8) \ 83066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 9) \ 84066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 10) 85066ea43fSLisandro Dalcin 86066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(VTX, 0) 87066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(SEG) 88066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TRI) 89066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(QUA) 90066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TET) 91066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(HEX) 92066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PRI) 93066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PYR) 94066ea43fSLisandro Dalcin 95066ea43fSLisandro Dalcin typedef enum { 96066ea43fSLisandro Dalcin GMSH_VTX = 0, 97066ea43fSLisandro Dalcin GMSH_SEG = 1, 98066ea43fSLisandro Dalcin GMSH_TRI = 2, 99066ea43fSLisandro Dalcin GMSH_QUA = 3, 100066ea43fSLisandro Dalcin GMSH_TET = 4, 101066ea43fSLisandro Dalcin GMSH_HEX = 5, 102066ea43fSLisandro Dalcin GMSH_PRI = 6, 103066ea43fSLisandro Dalcin GMSH_PYR = 7, 104066ea43fSLisandro Dalcin GMSH_NUM_POLYTOPES = 8 105066ea43fSLisandro Dalcin } GmshPolytopeType; 106066ea43fSLisandro Dalcin 107de78e4feSLisandro Dalcin typedef struct { 1080598e1a2SLisandro Dalcin int cellType; 109066ea43fSLisandro Dalcin int polytope; 1100598e1a2SLisandro Dalcin int dim; 1110598e1a2SLisandro Dalcin int order; 112066ea43fSLisandro Dalcin int numVerts; 1130598e1a2SLisandro Dalcin int numNodes; 114066ea43fSLisandro Dalcin int *(*lexorder)(void); 1150598e1a2SLisandro Dalcin } GmshCellInfo; 1160598e1a2SLisandro Dalcin 117066ea43fSLisandro Dalcin #define GmshCellEntry(cellType, polytope, dim, order) \ 1189371c9d4SSatish Balay { cellType, GMSH_##polytope, dim, order, GmshNumNodes_##polytope(1), GmshNumNodes_##polytope(order), GmshLexOrder_##polytope##_##order } 1190598e1a2SLisandro Dalcin 1200598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = { 121066ea43fSLisandro Dalcin GmshCellEntry(15, VTX, 0, 0), 1220598e1a2SLisandro Dalcin 1239371c9d4SSatish Balay GmshCellEntry(1, SEG, 1, 1), GmshCellEntry(8, SEG, 1, 2), GmshCellEntry(26, SEG, 1, 3), 1249371c9d4SSatish Balay GmshCellEntry(27, SEG, 1, 4), GmshCellEntry(28, SEG, 1, 5), GmshCellEntry(62, SEG, 1, 6), 1259371c9d4SSatish Balay GmshCellEntry(63, SEG, 1, 7), GmshCellEntry(64, SEG, 1, 8), GmshCellEntry(65, SEG, 1, 9), 126066ea43fSLisandro Dalcin GmshCellEntry(66, SEG, 1, 10), 1270598e1a2SLisandro Dalcin 1289371c9d4SSatish Balay GmshCellEntry(2, TRI, 2, 1), GmshCellEntry(9, TRI, 2, 2), GmshCellEntry(21, TRI, 2, 3), 1299371c9d4SSatish Balay GmshCellEntry(23, TRI, 2, 4), GmshCellEntry(25, TRI, 2, 5), GmshCellEntry(42, TRI, 2, 6), 1309371c9d4SSatish Balay GmshCellEntry(43, TRI, 2, 7), GmshCellEntry(44, TRI, 2, 8), GmshCellEntry(45, TRI, 2, 9), 131066ea43fSLisandro Dalcin GmshCellEntry(46, TRI, 2, 10), 1320598e1a2SLisandro Dalcin 1339371c9d4SSatish Balay GmshCellEntry(3, QUA, 2, 1), GmshCellEntry(10, QUA, 2, 2), {16, GMSH_QUA, 2, 2, 4, 8, GmshLexOrder_QUA_2_Serendipity}, 1349371c9d4SSatish Balay GmshCellEntry(36, QUA, 2, 3), GmshCellEntry(37, QUA, 2, 4), GmshCellEntry(38, QUA, 2, 5), 1359371c9d4SSatish Balay GmshCellEntry(47, QUA, 2, 6), GmshCellEntry(48, QUA, 2, 7), GmshCellEntry(49, QUA, 2, 8), 1369371c9d4SSatish Balay GmshCellEntry(50, QUA, 2, 9), GmshCellEntry(51, QUA, 2, 10), 1370598e1a2SLisandro Dalcin 1389371c9d4SSatish Balay GmshCellEntry(4, TET, 3, 1), GmshCellEntry(11, TET, 3, 2), GmshCellEntry(29, TET, 3, 3), 1399371c9d4SSatish Balay GmshCellEntry(30, TET, 3, 4), GmshCellEntry(31, TET, 3, 5), GmshCellEntry(71, TET, 3, 6), 1409371c9d4SSatish Balay GmshCellEntry(72, TET, 3, 7), GmshCellEntry(73, TET, 3, 8), GmshCellEntry(74, TET, 3, 9), 141066ea43fSLisandro Dalcin GmshCellEntry(75, TET, 3, 10), 1420598e1a2SLisandro Dalcin 1439371c9d4SSatish Balay GmshCellEntry(5, HEX, 3, 1), GmshCellEntry(12, HEX, 3, 2), {17, GMSH_HEX, 3, 2, 8, 20, GmshLexOrder_HEX_2_Serendipity}, 1449371c9d4SSatish Balay GmshCellEntry(92, HEX, 3, 3), GmshCellEntry(93, HEX, 3, 4), GmshCellEntry(94, HEX, 3, 5), 1459371c9d4SSatish Balay GmshCellEntry(95, HEX, 3, 6), GmshCellEntry(96, HEX, 3, 7), GmshCellEntry(97, HEX, 3, 8), 1469371c9d4SSatish Balay GmshCellEntry(98, HEX, 3, 9), GmshCellEntry(-1, HEX, 3, 10), 1470598e1a2SLisandro Dalcin 1489371c9d4SSatish Balay GmshCellEntry(6, PRI, 3, 1), GmshCellEntry(13, PRI, 3, 2), GmshCellEntry(90, PRI, 3, 3), 1499371c9d4SSatish Balay GmshCellEntry(91, PRI, 3, 4), GmshCellEntry(106, PRI, 3, 5), GmshCellEntry(107, PRI, 3, 6), 1509371c9d4SSatish Balay GmshCellEntry(108, PRI, 3, 7), GmshCellEntry(109, PRI, 3, 8), GmshCellEntry(110, PRI, 3, 9), 151066ea43fSLisandro Dalcin GmshCellEntry(-1, PRI, 3, 10), 1520598e1a2SLisandro Dalcin 1539371c9d4SSatish Balay GmshCellEntry(7, PYR, 3, 1), GmshCellEntry(14, PYR, 3, 2), GmshCellEntry(118, PYR, 3, 3), 1549371c9d4SSatish Balay GmshCellEntry(119, PYR, 3, 4), GmshCellEntry(120, PYR, 3, 5), GmshCellEntry(121, PYR, 3, 6), 1559371c9d4SSatish Balay GmshCellEntry(122, PYR, 3, 7), GmshCellEntry(123, PYR, 3, 8), GmshCellEntry(124, PYR, 3, 9), 156066ea43fSLisandro Dalcin GmshCellEntry(-1, PYR, 3, 10) 1570598e1a2SLisandro Dalcin 1580598e1a2SLisandro Dalcin #if 0 159066ea43fSLisandro Dalcin {20, GMSH_TRI, 2, 3, 3, 9, NULL}, 160066ea43fSLisandro Dalcin {18, GMSH_PRI, 3, 2, 6, 15, NULL}, 161066ea43fSLisandro Dalcin {19, GMSH_PYR, 3, 2, 5, 13, NULL}, 1620598e1a2SLisandro Dalcin #endif 1630598e1a2SLisandro Dalcin }; 1640598e1a2SLisandro Dalcin 1650598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150]; 1660598e1a2SLisandro Dalcin 1679371c9d4SSatish Balay static PetscErrorCode GmshCellInfoSetUp(void) { 1680598e1a2SLisandro Dalcin size_t i, n; 1690598e1a2SLisandro Dalcin static PetscBool called = PETSC_FALSE; 1700598e1a2SLisandro Dalcin 1710598e1a2SLisandro Dalcin if (called) return 0; 1720598e1a2SLisandro Dalcin PetscFunctionBegin; 1730598e1a2SLisandro Dalcin called = PETSC_TRUE; 174dd39110bSPierre Jolivet n = PETSC_STATIC_ARRAY_LENGTH(GmshCellMap); 1750598e1a2SLisandro Dalcin for (i = 0; i < n; ++i) { 1760598e1a2SLisandro Dalcin GmshCellMap[i].cellType = -1; 177066ea43fSLisandro Dalcin GmshCellMap[i].polytope = -1; 1780598e1a2SLisandro Dalcin } 179dd39110bSPierre Jolivet n = PETSC_STATIC_ARRAY_LENGTH(GmshCellTable); 180066ea43fSLisandro Dalcin for (i = 0; i < n; ++i) { 181066ea43fSLisandro Dalcin if (GmshCellTable[i].cellType <= 0) continue; 182066ea43fSLisandro Dalcin GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i]; 183066ea43fSLisandro Dalcin } 1840598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1850598e1a2SLisandro Dalcin } 1860598e1a2SLisandro Dalcin 1879371c9d4SSatish Balay #define GmshCellTypeCheck(ct) \ 1889371c9d4SSatish Balay 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_); \ 1899371c9d4SSatish Balay PetscCheck(GmshCellMap[_ct_].polytope != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_);) 1900598e1a2SLisandro Dalcin 1910598e1a2SLisandro Dalcin typedef struct { 192698a79b9SLisandro Dalcin PetscViewer viewer; 193698a79b9SLisandro Dalcin int fileFormat; 194698a79b9SLisandro Dalcin int dataSize; 195698a79b9SLisandro Dalcin PetscBool binary; 196698a79b9SLisandro Dalcin PetscBool byteSwap; 197698a79b9SLisandro Dalcin size_t wlen; 198698a79b9SLisandro Dalcin void *wbuf; 199698a79b9SLisandro Dalcin size_t slen; 200698a79b9SLisandro Dalcin void *sbuf; 2010598e1a2SLisandro Dalcin PetscInt *nbuf; 2020598e1a2SLisandro Dalcin PetscInt nodeStart; 2030598e1a2SLisandro Dalcin PetscInt nodeEnd; 20433a088b6SMatthew G. Knepley PetscInt *nodeMap; 205698a79b9SLisandro Dalcin } GmshFile; 206de78e4feSLisandro Dalcin 2079371c9d4SSatish Balay static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) { 208de78e4feSLisandro Dalcin size_t size = count * eltsize; 20911c56cb0SLisandro Dalcin 21011c56cb0SLisandro Dalcin PetscFunctionBegin; 211698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 2129566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->wbuf)); 2139566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size, &gmsh->wbuf)); 214698a79b9SLisandro Dalcin gmsh->wlen = size; 215de78e4feSLisandro Dalcin } 216698a79b9SLisandro Dalcin *(void **)buf = size ? gmsh->wbuf : NULL; 217de78e4feSLisandro Dalcin PetscFunctionReturn(0); 218de78e4feSLisandro Dalcin } 219de78e4feSLisandro Dalcin 2209371c9d4SSatish Balay static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) { 221698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 222698a79b9SLisandro Dalcin size_t size = count * dataSize; 223698a79b9SLisandro Dalcin 224698a79b9SLisandro Dalcin PetscFunctionBegin; 225698a79b9SLisandro Dalcin if (gmsh->slen < size) { 2269566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->sbuf)); 2279566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size, &gmsh->sbuf)); 228698a79b9SLisandro Dalcin gmsh->slen = size; 229698a79b9SLisandro Dalcin } 230698a79b9SLisandro Dalcin *(void **)buf = size ? gmsh->sbuf : NULL; 231698a79b9SLisandro Dalcin PetscFunctionReturn(0); 232698a79b9SLisandro Dalcin } 233698a79b9SLisandro Dalcin 2349371c9d4SSatish Balay static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) { 235de78e4feSLisandro Dalcin PetscFunctionBegin; 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype)); 2379566063dSJacob Faibussowitsch if (gmsh->byteSwap) PetscCall(PetscByteSwap(buf, dtype, count)); 238698a79b9SLisandro Dalcin PetscFunctionReturn(0); 239698a79b9SLisandro Dalcin } 240698a79b9SLisandro Dalcin 2419371c9d4SSatish Balay static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) { 242698a79b9SLisandro Dalcin PetscFunctionBegin; 2439566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING)); 244698a79b9SLisandro Dalcin PetscFunctionReturn(0); 245698a79b9SLisandro Dalcin } 246698a79b9SLisandro Dalcin 2479371c9d4SSatish Balay static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) { 248d6689ca9SLisandro Dalcin PetscFunctionBegin; 2499566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(line, Section, match)); 250d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 251d6689ca9SLisandro Dalcin } 252d6689ca9SLisandro Dalcin 2539371c9d4SSatish Balay static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) { 254d6689ca9SLisandro Dalcin PetscBool match; 255d6689ca9SLisandro Dalcin 256d6689ca9SLisandro Dalcin PetscFunctionBegin; 2579566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, Section, line, &match)); 25828b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s", Section); 259d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 260d6689ca9SLisandro Dalcin } 261d6689ca9SLisandro Dalcin 2629371c9d4SSatish Balay static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) { 263d6689ca9SLisandro Dalcin PetscBool match; 264d6689ca9SLisandro Dalcin 265d6689ca9SLisandro Dalcin PetscFunctionBegin; 266d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 2679566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 2689566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$Comments", line, &match)); 269d6689ca9SLisandro Dalcin if (!match) break; 270d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 2719566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 2729566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$EndComments", line, &match)); 273d6689ca9SLisandro Dalcin if (match) break; 274d6689ca9SLisandro Dalcin } 275d6689ca9SLisandro Dalcin } 276d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 277d6689ca9SLisandro Dalcin } 278d6689ca9SLisandro Dalcin 2799371c9d4SSatish Balay static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) { 280d6689ca9SLisandro Dalcin PetscFunctionBegin; 2819566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 2829566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, EndSection, line)); 283d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 284d6689ca9SLisandro Dalcin } 285d6689ca9SLisandro Dalcin 2869371c9d4SSatish Balay static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) { 287698a79b9SLisandro Dalcin PetscInt i; 288698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 289698a79b9SLisandro Dalcin 290698a79b9SLisandro Dalcin PetscFunctionBegin; 291698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 2929566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_INT)); 293698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 294698a79b9SLisandro Dalcin int *ibuf = NULL; 2959566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 2969566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_ENUM)); 297698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 298698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 299698a79b9SLisandro Dalcin long *ibuf = NULL; 3009566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 3019566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_LONG)); 302698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 303698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 304698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 3059566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 3069566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_INT64)); 307698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 308698a79b9SLisandro Dalcin } 309698a79b9SLisandro Dalcin PetscFunctionReturn(0); 310698a79b9SLisandro Dalcin } 311698a79b9SLisandro Dalcin 3129371c9d4SSatish Balay static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) { 313698a79b9SLisandro Dalcin PetscFunctionBegin; 3149566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_ENUM)); 315698a79b9SLisandro Dalcin PetscFunctionReturn(0); 316698a79b9SLisandro Dalcin } 317698a79b9SLisandro Dalcin 3189371c9d4SSatish Balay static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) { 319698a79b9SLisandro Dalcin PetscFunctionBegin; 3209566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_DOUBLE)); 321de78e4feSLisandro Dalcin PetscFunctionReturn(0); 322de78e4feSLisandro Dalcin } 323de78e4feSLisandro Dalcin 3249c0edc59SMatthew G. Knepley #define GMSH_MAX_TAGS 16 3257d5b9244SMatthew G. Knepley 326de78e4feSLisandro Dalcin typedef struct { 3270598e1a2SLisandro Dalcin PetscInt id; /* Entity ID */ 3280598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 329de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 330de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 3317d5b9244SMatthew G. Knepley int tags[GMSH_MAX_TAGS]; /* Tag array */ 332de78e4feSLisandro Dalcin } GmshEntity; 333de78e4feSLisandro Dalcin 334de78e4feSLisandro Dalcin typedef struct { 335730356f6SLisandro Dalcin GmshEntity *entity[4]; 336730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 337730356f6SLisandro Dalcin } GmshEntities; 338730356f6SLisandro Dalcin 3399371c9d4SSatish Balay static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) { 340730356f6SLisandro Dalcin PetscInt dim; 341730356f6SLisandro Dalcin 342730356f6SLisandro Dalcin PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(PetscNew(entities)); 344730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 3459566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(count[dim], &(*entities)->entity[dim])); 3469566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&(*entities)->entityMap[dim])); 347730356f6SLisandro Dalcin } 348730356f6SLisandro Dalcin PetscFunctionReturn(0); 349730356f6SLisandro Dalcin } 350730356f6SLisandro Dalcin 3519371c9d4SSatish Balay static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) { 3520598e1a2SLisandro Dalcin PetscInt dim; 3530598e1a2SLisandro Dalcin 3540598e1a2SLisandro Dalcin PetscFunctionBegin; 3550598e1a2SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 3560598e1a2SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 3579566063dSJacob Faibussowitsch PetscCall(PetscFree((*entities)->entity[dim])); 3589566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&(*entities)->entityMap[dim])); 3590598e1a2SLisandro Dalcin } 3609566063dSJacob Faibussowitsch PetscCall(PetscFree((*entities))); 3610598e1a2SLisandro Dalcin PetscFunctionReturn(0); 3620598e1a2SLisandro Dalcin } 3630598e1a2SLisandro Dalcin 3649371c9d4SSatish Balay static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity **entity) { 365730356f6SLisandro Dalcin PetscFunctionBegin; 3669566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(entities->entityMap[dim], eid, index)); 367730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 368730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 369730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 370730356f6SLisandro Dalcin PetscFunctionReturn(0); 371730356f6SLisandro Dalcin } 372730356f6SLisandro Dalcin 3739371c9d4SSatish Balay static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity **entity) { 374730356f6SLisandro Dalcin PetscInt index; 375730356f6SLisandro Dalcin 376730356f6SLisandro Dalcin PetscFunctionBegin; 3779566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(entities->entityMap[dim], eid, &index)); 378730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 379730356f6SLisandro Dalcin PetscFunctionReturn(0); 380730356f6SLisandro Dalcin } 381730356f6SLisandro Dalcin 3820598e1a2SLisandro Dalcin typedef struct { 3830598e1a2SLisandro Dalcin PetscInt *id; /* Node IDs */ 3840598e1a2SLisandro Dalcin double *xyz; /* Coordinates */ 38581a1af93SMatthew G. Knepley PetscInt *tag; /* Physical tag */ 3860598e1a2SLisandro Dalcin } GmshNodes; 3870598e1a2SLisandro Dalcin 3889371c9d4SSatish Balay static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes) { 389730356f6SLisandro Dalcin PetscFunctionBegin; 3909566063dSJacob Faibussowitsch PetscCall(PetscNew(nodes)); 3919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count * 1, &(*nodes)->id)); 3929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count * 3, &(*nodes)->xyz)); 3937d5b9244SMatthew G. Knepley PetscCall(PetscMalloc1(count * GMSH_MAX_TAGS, &(*nodes)->tag)); 3940598e1a2SLisandro Dalcin PetscFunctionReturn(0); 395730356f6SLisandro Dalcin } 3960598e1a2SLisandro Dalcin 3979371c9d4SSatish Balay static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes) { 3980598e1a2SLisandro Dalcin PetscFunctionBegin; 3990598e1a2SLisandro Dalcin if (!*nodes) PetscFunctionReturn(0); 4009566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->id)); 4019566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->xyz)); 4029566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->tag)); 4039566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes))); 404730356f6SLisandro Dalcin PetscFunctionReturn(0); 405730356f6SLisandro Dalcin } 406730356f6SLisandro Dalcin 407730356f6SLisandro Dalcin typedef struct { 4080598e1a2SLisandro Dalcin PetscInt id; /* Element ID */ 4090598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 410de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 4110598e1a2SLisandro Dalcin PetscInt numVerts; /* Size of vertex array */ 412de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 4130598e1a2SLisandro Dalcin PetscInt *nodes; /* Vertex/Node array */ 41481a1af93SMatthew G. Knepley PetscInt numTags; /* Size of physical tag array */ 4157d5b9244SMatthew G. Knepley int tags[GMSH_MAX_TAGS]; /* Physical tag array */ 416de78e4feSLisandro Dalcin } GmshElement; 417de78e4feSLisandro Dalcin 4189371c9d4SSatish Balay static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements) { 4190598e1a2SLisandro Dalcin PetscFunctionBegin; 4209566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(count, elements)); 4210598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4220598e1a2SLisandro Dalcin } 4230598e1a2SLisandro Dalcin 4249371c9d4SSatish Balay static PetscErrorCode GmshElementsDestroy(GmshElement **elements) { 4250598e1a2SLisandro Dalcin PetscFunctionBegin; 4260598e1a2SLisandro Dalcin if (!*elements) PetscFunctionReturn(0); 4279566063dSJacob Faibussowitsch PetscCall(PetscFree(*elements)); 4280598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4290598e1a2SLisandro Dalcin } 4300598e1a2SLisandro Dalcin 4310598e1a2SLisandro Dalcin typedef struct { 4320598e1a2SLisandro Dalcin PetscInt dim; 433066ea43fSLisandro Dalcin PetscInt order; 4340598e1a2SLisandro Dalcin GmshEntities *entities; 4350598e1a2SLisandro Dalcin PetscInt numNodes; 4360598e1a2SLisandro Dalcin GmshNodes *nodelist; 4370598e1a2SLisandro Dalcin PetscInt numElems; 4380598e1a2SLisandro Dalcin GmshElement *elements; 4390598e1a2SLisandro Dalcin PetscInt numVerts; 4400598e1a2SLisandro Dalcin PetscInt numCells; 4410598e1a2SLisandro Dalcin PetscInt *periodMap; 4420598e1a2SLisandro Dalcin PetscInt *vertexMap; 4430598e1a2SLisandro Dalcin PetscSegBuffer segbuf; 444a45dabc8SMatthew G. Knepley PetscInt numRegions; 44590d778b8SLisandro Dalcin PetscInt *regionDims; 446a45dabc8SMatthew G. Knepley PetscInt *regionTags; 447a45dabc8SMatthew G. Knepley char **regionNames; 4480598e1a2SLisandro Dalcin } GmshMesh; 4490598e1a2SLisandro Dalcin 4509371c9d4SSatish Balay static PetscErrorCode GmshMeshCreate(GmshMesh **mesh) { 4510598e1a2SLisandro Dalcin PetscFunctionBegin; 4529566063dSJacob Faibussowitsch PetscCall(PetscNew(mesh)); 4539566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf)); 4540598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4550598e1a2SLisandro Dalcin } 4560598e1a2SLisandro Dalcin 4579371c9d4SSatish Balay static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh) { 458a45dabc8SMatthew G. Knepley PetscInt r; 4590598e1a2SLisandro Dalcin 4600598e1a2SLisandro Dalcin PetscFunctionBegin; 4610598e1a2SLisandro Dalcin if (!*mesh) PetscFunctionReturn(0); 4629566063dSJacob Faibussowitsch PetscCall(GmshEntitiesDestroy(&(*mesh)->entities)); 4639566063dSJacob Faibussowitsch PetscCall(GmshNodesDestroy(&(*mesh)->nodelist)); 4649566063dSJacob Faibussowitsch PetscCall(GmshElementsDestroy(&(*mesh)->elements)); 4659566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh)->periodMap)); 4669566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh)->vertexMap)); 4679566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&(*mesh)->segbuf)); 4689566063dSJacob Faibussowitsch for (r = 0; r < (*mesh)->numRegions; ++r) PetscCall(PetscFree((*mesh)->regionNames[r])); 46990d778b8SLisandro Dalcin PetscCall(PetscFree3((*mesh)->regionDims, (*mesh)->regionTags, (*mesh)->regionNames)); 4709566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh))); 4710598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4720598e1a2SLisandro Dalcin } 4730598e1a2SLisandro Dalcin 4749371c9d4SSatish Balay static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh) { 475698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 476698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 477de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 4787d5b9244SMatthew G. Knepley int n, t, num, nid, snum; 4790598e1a2SLisandro Dalcin GmshNodes *nodes; 480de78e4feSLisandro Dalcin 481de78e4feSLisandro Dalcin PetscFunctionBegin; 4829566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 483698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 48408401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 4859566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(num, &nodes)); 4860598e1a2SLisandro Dalcin mesh->numNodes = num; 4870598e1a2SLisandro Dalcin mesh->nodelist = nodes; 4880598e1a2SLisandro Dalcin for (n = 0; n < num; ++n) { 4890598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n * 3; 4909566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM)); 4919566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE)); 4929566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 4939566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 4940598e1a2SLisandro Dalcin nodes->id[n] = nid; 4957d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1; 49611c56cb0SLisandro Dalcin } 49711c56cb0SLisandro Dalcin PetscFunctionReturn(0); 49811c56cb0SLisandro Dalcin } 49911c56cb0SLisandro Dalcin 500de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 501de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 502de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 503de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 5049371c9d4SSatish Balay static PetscErrorCode GmshReadElements_v22(GmshFile *gmsh, GmshMesh *mesh) { 505698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 506698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 507698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 508de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 5090598e1a2SLisandro Dalcin int i, c, p, num, ibuf[1 + 4 + 1000], snum; 5100598e1a2SLisandro Dalcin int cellType, numElem, numVerts, numNodes, numTags; 511df032b82SMatthew G. Knepley GmshElement *elements; 5120598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 513df032b82SMatthew G. Knepley 514df032b82SMatthew G. Knepley PetscFunctionBegin; 5159566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 516698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 51708401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 5189566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(num, &elements)); 5190598e1a2SLisandro Dalcin mesh->numElems = num; 5200598e1a2SLisandro Dalcin mesh->elements = elements; 521698a79b9SLisandro Dalcin for (c = 0; c < num;) { 5229566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM)); 5239566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3)); 5240598e1a2SLisandro Dalcin 5250598e1a2SLisandro Dalcin cellType = binary ? ibuf[0] : ibuf[1]; 5260598e1a2SLisandro Dalcin numElem = binary ? ibuf[1] : 1; 527df032b82SMatthew G. Knepley numTags = ibuf[2]; 5280598e1a2SLisandro Dalcin 5299566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 5300598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 5310598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 5320598e1a2SLisandro Dalcin 533df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 5340598e1a2SLisandro Dalcin GmshElement *element = elements + c; 5350598e1a2SLisandro Dalcin const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off; 5360598e1a2SLisandro Dalcin const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1; 5379566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf + off, nint, NULL, PETSC_ENUM)); 5389566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf + off, PETSC_ENUM, nint)); 5390598e1a2SLisandro Dalcin element->id = id[0]; 5400598e1a2SLisandro Dalcin element->dim = GmshCellMap[cellType].dim; 5410598e1a2SLisandro Dalcin element->cellType = cellType; 5420598e1a2SLisandro Dalcin element->numVerts = numVerts; 5430598e1a2SLisandro Dalcin element->numNodes = numNodes; 5447d5b9244SMatthew G. Knepley element->numTags = PetscMin(numTags, GMSH_MAX_TAGS); 5459566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 5460598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 5470598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = tags[p]; 548df032b82SMatthew G. Knepley } 549df032b82SMatthew G. Knepley } 550df032b82SMatthew G. Knepley PetscFunctionReturn(0); 551df032b82SMatthew G. Knepley } 5527d282ae0SMichael Lange 553de78e4feSLisandro Dalcin /* 554de78e4feSLisandro Dalcin $Entities 555de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 556de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 557de78e4feSLisandro Dalcin // points 558de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 559de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 560de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 561de78e4feSLisandro Dalcin ... 562de78e4feSLisandro Dalcin // curves 563de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 564de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 565de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 566de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 567de78e4feSLisandro Dalcin ... 568de78e4feSLisandro Dalcin // surfaces 569de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 570de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 571de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 572de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 573de78e4feSLisandro Dalcin ... 574de78e4feSLisandro Dalcin // volumes 575de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 576de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 577de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 578de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 579de78e4feSLisandro Dalcin ... 580de78e4feSLisandro Dalcin $EndEntities 581de78e4feSLisandro Dalcin */ 5829371c9d4SSatish Balay static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh) { 583698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 584698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 585698a79b9SLisandro Dalcin long index, num, lbuf[4]; 586730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 587698a79b9SLisandro Dalcin PetscInt count[4], i; 588a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 589de78e4feSLisandro Dalcin 590de78e4feSLisandro Dalcin PetscFunctionBegin; 5919566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG)); 5929566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(lbuf, PETSC_LONG, 4)); 593698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 5949566063dSJacob Faibussowitsch PetscCall(GmshEntitiesCreate(count, &mesh->entities)); 595de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 596730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 5979566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM)); 5989566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&eid, PETSC_ENUM, 1)); 5999566063dSJacob Faibussowitsch PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity)); 6009566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE)); 6019566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6)); 6029566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG)); 6039566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1)); 6049566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf)); 6059566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM)); 6069566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num)); 6077d5b9244SMatthew G. Knepley entity->numTags = numTags = (int)PetscMin(num, GMSH_MAX_TAGS); 608de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 609de78e4feSLisandro Dalcin if (dim == 0) continue; 6109566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG)); 6119566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1)); 6129566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf)); 6139566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM)); 6149566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num)); 615de78e4feSLisandro Dalcin } 616de78e4feSLisandro Dalcin } 617de78e4feSLisandro Dalcin PetscFunctionReturn(0); 618de78e4feSLisandro Dalcin } 619de78e4feSLisandro Dalcin 620de78e4feSLisandro Dalcin /* 621de78e4feSLisandro Dalcin $Nodes 622de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 623de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 624de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 625de78e4feSLisandro Dalcin ... 626de78e4feSLisandro Dalcin ... 627de78e4feSLisandro Dalcin $EndNodes 628de78e4feSLisandro Dalcin */ 6299371c9d4SSatish Balay static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh) { 630698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 631698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 6327d5b9244SMatthew G. Knepley long block, node, n, t, numEntityBlocks, numTotalNodes, numNodes; 633de78e4feSLisandro Dalcin int info[3], nid; 6340598e1a2SLisandro Dalcin GmshNodes *nodes; 635de78e4feSLisandro Dalcin 636de78e4feSLisandro Dalcin PetscFunctionBegin; 6379566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG)); 6389566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1)); 6399566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG)); 6409566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numTotalNodes, PETSC_LONG, 1)); 6419566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(numTotalNodes, &nodes)); 6420598e1a2SLisandro Dalcin mesh->numNodes = numTotalNodes; 6430598e1a2SLisandro Dalcin mesh->nodelist = nodes; 6440598e1a2SLisandro Dalcin for (n = 0, block = 0; block < numEntityBlocks; ++block) { 6459566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM)); 6469566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG)); 6479566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numNodes, PETSC_LONG, 1)); 648698a79b9SLisandro Dalcin if (gmsh->binary) { 649277f51e8SBarry Smith size_t nbytes = sizeof(int) + 3 * sizeof(double); 650277f51e8SBarry Smith char *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */ 6519566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numNodes, nbytes, &cbuf)); 6529566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, cbuf, numNodes * nbytes, NULL, PETSC_CHAR)); 6530598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 654de78e4feSLisandro Dalcin char *cnid = cbuf + node * nbytes, *cxyz = cnid + sizeof(int); 6550598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n * 3; 6569566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cnid, PETSC_ENUM, 1)); 6579566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cxyz, PETSC_DOUBLE, 3)); 6589566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&nid, cnid, sizeof(int))); 6599566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(xyz, cxyz, 3 * sizeof(double))); 6609566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 6619566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 6620598e1a2SLisandro Dalcin nodes->id[n] = nid; 6637d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1; 664de78e4feSLisandro Dalcin } 665de78e4feSLisandro Dalcin } else { 6660598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 6670598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n * 3; 6689566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM)); 6699566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE)); 6709566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 6719566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 6720598e1a2SLisandro Dalcin nodes->id[n] = nid; 6737d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n * GMSH_MAX_TAGS + t] = -1; 674de78e4feSLisandro Dalcin } 675de78e4feSLisandro Dalcin } 676de78e4feSLisandro Dalcin } 677de78e4feSLisandro Dalcin PetscFunctionReturn(0); 678de78e4feSLisandro Dalcin } 679de78e4feSLisandro Dalcin 680de78e4feSLisandro Dalcin /* 681de78e4feSLisandro Dalcin $Elements 682de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 683de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 684de78e4feSLisandro Dalcin tag(int) numVert[...](int) 685de78e4feSLisandro Dalcin ... 686de78e4feSLisandro Dalcin ... 687de78e4feSLisandro Dalcin $EndElements 688de78e4feSLisandro Dalcin */ 6899371c9d4SSatish Balay static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh) { 690698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 691698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 692de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 693de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 6940598e1a2SLisandro Dalcin int eid, dim, cellType, numVerts, numNodes, numTags; 695a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 696de78e4feSLisandro Dalcin GmshElement *elements; 6970598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 698de78e4feSLisandro Dalcin 699de78e4feSLisandro Dalcin PetscFunctionBegin; 7009566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG)); 7019566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1)); 7029566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG)); 7039566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numTotalElements, PETSC_LONG, 1)); 7049566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(numTotalElements, &elements)); 7050598e1a2SLisandro Dalcin mesh->numElems = numTotalElements; 7060598e1a2SLisandro Dalcin mesh->elements = elements; 707de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 7089566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM)); 7099566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(info, PETSC_ENUM, 3)); 7109371c9d4SSatish Balay eid = info[0]; 7119371c9d4SSatish Balay dim = info[1]; 7129371c9d4SSatish Balay cellType = info[2]; 7139566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 7149566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 7150598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 7160598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 717730356f6SLisandro Dalcin numTags = entity->numTags; 7189566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG)); 7199566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numElements, PETSC_LONG, 1)); 7209566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, (1 + numNodes) * numElements, sizeof(int), &ibuf)); 7219566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, (1 + numNodes) * numElements, NULL, PETSC_ENUM)); 7229566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, (1 + numNodes) * numElements)); 723de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 724de78e4feSLisandro Dalcin GmshElement *element = elements + c; 7250598e1a2SLisandro Dalcin const int *id = ibuf + elem * (1 + numNodes), *nodes = id + 1; 7260598e1a2SLisandro Dalcin element->id = id[0]; 727de78e4feSLisandro Dalcin element->dim = dim; 728de78e4feSLisandro Dalcin element->cellType = cellType; 7290598e1a2SLisandro Dalcin element->numVerts = numVerts; 730de78e4feSLisandro Dalcin element->numNodes = numNodes; 731de78e4feSLisandro Dalcin element->numTags = numTags; 7329566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 7330598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 7340598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 735de78e4feSLisandro Dalcin } 736de78e4feSLisandro Dalcin } 737698a79b9SLisandro Dalcin PetscFunctionReturn(0); 738698a79b9SLisandro Dalcin } 739698a79b9SLisandro Dalcin 7409371c9d4SSatish Balay static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[]) { 741698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 742698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 743698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 744698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 745698a79b9SLisandro Dalcin int numPeriodic, snum, i; 746698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 7470598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 748698a79b9SLisandro Dalcin 749698a79b9SLisandro Dalcin PetscFunctionBegin; 750698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 7519566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 752698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 75308401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 754698a79b9SLisandro Dalcin } else { 7559566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM)); 7569566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numPeriodic, PETSC_ENUM, 1)); 757698a79b9SLisandro Dalcin } 758698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 7599dddd249SSatish Balay int ibuf[3], correspondingDim = -1, correspondingTag = -1, primaryTag = -1, correspondingNode, primaryNode; 760698a79b9SLisandro Dalcin long j, nNodes; 761698a79b9SLisandro Dalcin double affine[16]; 762698a79b9SLisandro Dalcin 763698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 7649566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING)); 7659dddd249SSatish Balay snum = sscanf(line, "%d %d %d", &correspondingDim, &correspondingTag, &primaryTag); 76608401ef6SPierre Jolivet PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 767698a79b9SLisandro Dalcin } else { 7689566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM)); 7699566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3)); 7709371c9d4SSatish Balay correspondingDim = ibuf[0]; 7719371c9d4SSatish Balay correspondingTag = ibuf[1]; 7729371c9d4SSatish Balay primaryTag = ibuf[2]; 773698a79b9SLisandro Dalcin } 7749371c9d4SSatish Balay (void)correspondingDim; 7759371c9d4SSatish Balay (void)correspondingTag; 7769371c9d4SSatish Balay (void)primaryTag; /* unused */ 777698a79b9SLisandro Dalcin 778698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 7799566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 780698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 7814c500f23SPierre Jolivet if (snum != 1) { /* discard transformation and try again */ 7829566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING)); 7839566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 784698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 78508401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 786698a79b9SLisandro Dalcin } 787698a79b9SLisandro Dalcin } else { 7889566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG)); 7899566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1)); 7904c500f23SPierre Jolivet if (nNodes == -1) { /* discard transformation and try again */ 7919566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE)); 7929566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG)); 7939566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1)); 794698a79b9SLisandro Dalcin } 795698a79b9SLisandro Dalcin } 796698a79b9SLisandro Dalcin 797698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 798698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 7999566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 8009dddd249SSatish Balay snum = sscanf(line, "%d %d", &correspondingNode, &primaryNode); 80108401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 802698a79b9SLisandro Dalcin } else { 8039566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM)); 8049566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 2)); 8059371c9d4SSatish Balay correspondingNode = ibuf[0]; 8069371c9d4SSatish Balay primaryNode = ibuf[1]; 807698a79b9SLisandro Dalcin } 8089dddd249SSatish Balay correspondingNode = (int)nodeMap[correspondingNode]; 8099dddd249SSatish Balay primaryNode = (int)nodeMap[primaryNode]; 8109dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 811698a79b9SLisandro Dalcin } 812698a79b9SLisandro Dalcin } 813698a79b9SLisandro Dalcin PetscFunctionReturn(0); 814698a79b9SLisandro Dalcin } 815698a79b9SLisandro Dalcin 8160598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 817698a79b9SLisandro Dalcin $Entities 818698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 819698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 820698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 821698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 822698a79b9SLisandro Dalcin ... 823698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 824698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 825698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 826698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 827698a79b9SLisandro Dalcin ... 828698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 829698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 830698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 831698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 832698a79b9SLisandro Dalcin ... 833698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 834698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 835698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 836698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 837698a79b9SLisandro Dalcin ... 838698a79b9SLisandro Dalcin $EndEntities 839698a79b9SLisandro Dalcin */ 8409371c9d4SSatish Balay static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh) { 841698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 842698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 843698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 844698a79b9SLisandro Dalcin 845698a79b9SLisandro Dalcin PetscFunctionBegin; 8469566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, count, 4)); 8479566063dSJacob Faibussowitsch PetscCall(GmshEntitiesCreate(count, &mesh->entities)); 848698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 849698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 8509566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, &eid, 1)); 8519566063dSJacob Faibussowitsch PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity)); 8529566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6)); 8539566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numTags, 1)); 8549566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags)); 8559566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, tags, numTags)); 8569c0edc59SMatthew G. Knepley PetscCheck(numTags <= GMSH_MAX_TAGS, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PETSc currently supports up to %" PetscInt_FMT " tags per entity, not %" PetscInt_FMT, (PetscInt)GMSH_MAX_TAGS, numTags); 8577d5b9244SMatthew G. Knepley entity->numTags = numTags; 858698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 859698a79b9SLisandro Dalcin if (dim == 0) continue; 8609566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numTags, 1)); 8619566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags)); 8629566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, tags, numTags)); 86381a1af93SMatthew G. Knepley /* Currently, we do not save the ids for the bounding entities */ 864698a79b9SLisandro Dalcin } 865698a79b9SLisandro Dalcin } 866698a79b9SLisandro Dalcin PetscFunctionReturn(0); 867698a79b9SLisandro Dalcin } 868698a79b9SLisandro Dalcin 86933a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 870698a79b9SLisandro Dalcin $Nodes 871698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 872698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 873698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 874698a79b9SLisandro Dalcin nodeTag(size_t) 875698a79b9SLisandro Dalcin ... 876698a79b9SLisandro Dalcin x(double) y(double) z(double) 877698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 878698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 879698a79b9SLisandro Dalcin ... 880698a79b9SLisandro Dalcin ... 881698a79b9SLisandro Dalcin $EndNodes 882698a79b9SLisandro Dalcin */ 8839371c9d4SSatish Balay static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh) { 88481a1af93SMatthew G. Knepley int info[3], dim, eid, parametric; 8857d5b9244SMatthew G. Knepley PetscInt sizes[4], numEntityBlocks, numTags, t, numNodes, numNodesBlock = 0, block, node, n; 88681a1af93SMatthew G. Knepley GmshEntity *entity = NULL; 8870598e1a2SLisandro Dalcin GmshNodes *nodes; 888698a79b9SLisandro Dalcin 889698a79b9SLisandro Dalcin PetscFunctionBegin; 8909566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, sizes, 4)); 8919371c9d4SSatish Balay numEntityBlocks = sizes[0]; 8929371c9d4SSatish Balay numNodes = sizes[1]; 8939566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(numNodes, &nodes)); 8940598e1a2SLisandro Dalcin mesh->numNodes = numNodes; 8950598e1a2SLisandro Dalcin mesh->nodelist = nodes; 896698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 8979566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 8989371c9d4SSatish Balay dim = info[0]; 8999371c9d4SSatish Balay eid = info[1]; 9009371c9d4SSatish Balay parametric = info[2]; 9019566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 9027d5b9244SMatthew G. Knepley numTags = entity->numTags; 90381a1af93SMatthew G. Knepley PetscCheck(!parametric, PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 9049566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numNodesBlock, 1)); 9059566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, nodes->id + node, numNodesBlock)); 9069566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, nodes->xyz + node * 3, numNodesBlock * 3)); 9077d5b9244SMatthew G. Knepley for (n = 0; n < numNodesBlock; ++n) { 9087d5b9244SMatthew G. Knepley PetscInt *tags = &nodes->tag[node * GMSH_MAX_TAGS]; 9097d5b9244SMatthew G. Knepley 9107d5b9244SMatthew G. Knepley for (t = 0; t < numTags; ++t) tags[n * GMSH_MAX_TAGS + t] = entity->tags[t]; 9117d5b9244SMatthew G. Knepley for (t = numTags; t < GMSH_MAX_TAGS; ++t) tags[n * GMSH_MAX_TAGS + t] = -1; 9127d5b9244SMatthew G. Knepley } 913698a79b9SLisandro Dalcin } 9140598e1a2SLisandro Dalcin gmsh->nodeStart = sizes[2]; 9150598e1a2SLisandro Dalcin gmsh->nodeEnd = sizes[3] + 1; 916698a79b9SLisandro Dalcin PetscFunctionReturn(0); 917698a79b9SLisandro Dalcin } 918698a79b9SLisandro Dalcin 91933a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 920698a79b9SLisandro Dalcin $Elements 921698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 922698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 923698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 924698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 925698a79b9SLisandro Dalcin ... 926698a79b9SLisandro Dalcin ... 927698a79b9SLisandro Dalcin $EndElements 928698a79b9SLisandro Dalcin */ 9299371c9d4SSatish Balay static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh) { 9300598e1a2SLisandro Dalcin int info[3], eid, dim, cellType; 9310598e1a2SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p; 932698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 933698a79b9SLisandro Dalcin GmshElement *elements; 9340598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 935698a79b9SLisandro Dalcin 936698a79b9SLisandro Dalcin PetscFunctionBegin; 9379566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, sizes, 4)); 9389371c9d4SSatish Balay numEntityBlocks = sizes[0]; 9399371c9d4SSatish Balay numElements = sizes[1]; 9409566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(numElements, &elements)); 9410598e1a2SLisandro Dalcin mesh->numElems = numElements; 9420598e1a2SLisandro Dalcin mesh->elements = elements; 943698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 9449566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 9459371c9d4SSatish Balay dim = info[0]; 9469371c9d4SSatish Balay eid = info[1]; 9479371c9d4SSatish Balay cellType = info[2]; 9489566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 9499566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 9500598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 9510598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 9527d5b9244SMatthew G. Knepley numTags = entity->numTags; 9539566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numBlockElements, 1)); 9549566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, (1 + numNodes) * numBlockElements, sizeof(PetscInt), &ibuf)); 9559566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, ibuf, (1 + numNodes) * numBlockElements)); 956698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 957698a79b9SLisandro Dalcin GmshElement *element = elements + c; 9580598e1a2SLisandro Dalcin const PetscInt *id = ibuf + elem * (1 + numNodes), *nodes = id + 1; 959698a79b9SLisandro Dalcin element->id = id[0]; 960698a79b9SLisandro Dalcin element->dim = dim; 961698a79b9SLisandro Dalcin element->cellType = cellType; 9620598e1a2SLisandro Dalcin element->numVerts = numVerts; 963698a79b9SLisandro Dalcin element->numNodes = numNodes; 964698a79b9SLisandro Dalcin element->numTags = numTags; 9659566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 9660598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 9670598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 968698a79b9SLisandro Dalcin } 969698a79b9SLisandro Dalcin } 970698a79b9SLisandro Dalcin PetscFunctionReturn(0); 971698a79b9SLisandro Dalcin } 972698a79b9SLisandro Dalcin 9730598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 974698a79b9SLisandro Dalcin $Periodic 975698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 9769dddd249SSatish Balay entityDim(int) entityTag(int) entityTagPrimary(int) 977698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 978698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 9799dddd249SSatish Balay nodeTag(size_t) nodeTagPrimary(size_t) 980698a79b9SLisandro Dalcin ... 981698a79b9SLisandro Dalcin ... 982698a79b9SLisandro Dalcin $EndPeriodic 983698a79b9SLisandro Dalcin */ 9849371c9d4SSatish Balay static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[]) { 985698a79b9SLisandro Dalcin int info[3]; 986698a79b9SLisandro Dalcin double dbuf[16]; 9870598e1a2SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 9880598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 989698a79b9SLisandro Dalcin 990698a79b9SLisandro Dalcin PetscFunctionBegin; 9919566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numPeriodicLinks, 1)); 992698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 9939566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 9949566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numAffine, 1)); 9959566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, dbuf, numAffine)); 9969566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numCorrespondingNodes, 1)); 9979566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags)); 9989566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, nodeTags, numCorrespondingNodes * 2)); 999698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 10009dddd249SSatish Balay PetscInt correspondingNode = nodeMap[nodeTags[node * 2 + 0]]; 10019dddd249SSatish Balay PetscInt primaryNode = nodeMap[nodeTags[node * 2 + 1]]; 10029dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 1003698a79b9SLisandro Dalcin } 1004698a79b9SLisandro Dalcin } 1005698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1006698a79b9SLisandro Dalcin } 1007698a79b9SLisandro Dalcin 10080598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 1009d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 1010d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 1011d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 1012d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 1013d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 1014d6689ca9SLisandro Dalcin $EndMeshFormat 1015d6689ca9SLisandro Dalcin */ 10169371c9d4SSatish Balay static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh) { 1017698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1018698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 1019698a79b9SLisandro Dalcin float version; 1020698a79b9SLisandro Dalcin 1021698a79b9SLisandro Dalcin PetscFunctionBegin; 10229566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 3)); 1023698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 1024a8d4e440SJunchao Zhang fileFormat = (int)roundf(version * 10); 102508401ef6SPierre Jolivet PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 1026a8d4e440SJunchao Zhang PetscCheck(fileFormat >= 22, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 10271dca8a05SBarry Smith PetscCheck((int)version != 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 1028a8d4e440SJunchao Zhang PetscCheck(fileFormat <= 41, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 102908401ef6SPierre Jolivet PetscCheck(!gmsh->binary || fileType, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII"); 103008401ef6SPierre Jolivet PetscCheck(gmsh->binary || !fileType, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary"); 10311dca8a05SBarry Smith PetscCheck(fileFormat > 40 || dataSize == sizeof(double), PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize); 10321dca8a05SBarry Smith 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); 1033698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 1034698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 1035698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 1036698a79b9SLisandro Dalcin if (gmsh->binary) { 10379566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, &checkEndian, 1)); 1038698a79b9SLisandro Dalcin if (checkEndian != 1) { 10399566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(&checkEndian, PETSC_ENUM, 1)); 104008401ef6SPierre Jolivet PetscCheck(checkEndian == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line); 1041698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 1042698a79b9SLisandro Dalcin } 1043698a79b9SLisandro Dalcin } 1044698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1045698a79b9SLisandro Dalcin } 1046698a79b9SLisandro Dalcin 10471f643af3SLisandro Dalcin /* 10481f643af3SLisandro Dalcin PhysicalNames 10491f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 10501f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 10511f643af3SLisandro Dalcin ... 10521f643af3SLisandro Dalcin $EndPhysicalNames 10531f643af3SLisandro Dalcin */ 10549371c9d4SSatish Balay static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh, GmshMesh *mesh) { 10555f5cd6d5SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN], name[128 + 2], *p, *q, *r; 1056a45dabc8SMatthew G. Knepley int snum, region, dim, tag; 1057698a79b9SLisandro Dalcin 1058698a79b9SLisandro Dalcin PetscFunctionBegin; 10599566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 1060a45dabc8SMatthew G. Knepley snum = sscanf(line, "%d", ®ion); 1061a45dabc8SMatthew G. Knepley mesh->numRegions = region; 106208401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 106390d778b8SLisandro Dalcin PetscCall(PetscMalloc3(mesh->numRegions, &mesh->regionDims, mesh->numRegions, &mesh->regionTags, mesh->numRegions, &mesh->regionNames)); 1064a45dabc8SMatthew G. Knepley for (region = 0; region < mesh->numRegions; ++region) { 10659566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 2)); 10661f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 106708401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10689566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, -(PetscInt)sizeof(line))); 10699566063dSJacob Faibussowitsch PetscCall(PetscStrchr(line, '"', &p)); 107028b400f6SJacob Faibussowitsch PetscCheck(p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10719566063dSJacob Faibussowitsch PetscCall(PetscStrrchr(line, '"', &q)); 107208401ef6SPierre Jolivet PetscCheck(q != p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10735f5cd6d5SMatthew G. Knepley PetscCall(PetscStrrchr(line, ':', &r)); 10745f5cd6d5SMatthew G. Knepley if (p != r) q = r; 10759566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, p + 1, (size_t)(q - p - 1))); 107690d778b8SLisandro Dalcin mesh->regionDims[region] = dim; 1077a45dabc8SMatthew G. Knepley mesh->regionTags[region] = tag; 10789566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &mesh->regionNames[region])); 1079698a79b9SLisandro Dalcin } 1080de78e4feSLisandro Dalcin PetscFunctionReturn(0); 1081de78e4feSLisandro Dalcin } 1082de78e4feSLisandro Dalcin 10839371c9d4SSatish Balay static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh) { 10840598e1a2SLisandro Dalcin PetscFunctionBegin; 10850598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 10869566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadEntities_v41(gmsh, mesh)); break; 10879566063dSJacob Faibussowitsch default: PetscCall(GmshReadEntities_v40(gmsh, mesh)); break; 108896ca5757SLisandro Dalcin } 10890598e1a2SLisandro Dalcin PetscFunctionReturn(0); 10900598e1a2SLisandro Dalcin } 10910598e1a2SLisandro Dalcin 10929371c9d4SSatish Balay static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh) { 10930598e1a2SLisandro Dalcin PetscFunctionBegin; 10940598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 10959566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadNodes_v41(gmsh, mesh)); break; 10969566063dSJacob Faibussowitsch case 40: PetscCall(GmshReadNodes_v40(gmsh, mesh)); break; 10979566063dSJacob Faibussowitsch default: PetscCall(GmshReadNodes_v22(gmsh, mesh)); break; 10980598e1a2SLisandro Dalcin } 10990598e1a2SLisandro Dalcin 11000598e1a2SLisandro Dalcin { /* Gmsh v2.2/v4.0 does not provide min/max node tags */ 11010598e1a2SLisandro Dalcin if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) { 11020598e1a2SLisandro Dalcin PetscInt tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n; 11030598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11040598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11050598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 11060598e1a2SLisandro Dalcin tagMin = PetscMin(tag, tagMin); 11070598e1a2SLisandro Dalcin tagMax = PetscMax(tag, tagMax); 11080598e1a2SLisandro Dalcin } 11090598e1a2SLisandro Dalcin gmsh->nodeStart = tagMin; 11100598e1a2SLisandro Dalcin gmsh->nodeEnd = tagMax + 1; 11110598e1a2SLisandro Dalcin } 11120598e1a2SLisandro Dalcin } 11130598e1a2SLisandro Dalcin 11140598e1a2SLisandro Dalcin { /* Support for sparse node tags */ 11150598e1a2SLisandro Dalcin PetscInt n, t; 11160598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf)); 11180598e1a2SLisandro Dalcin for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT; 11190598e1a2SLisandro Dalcin gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart; 11200598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11210598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 112263a3b9bcSJacob Faibussowitsch PetscCheck(gmsh->nodeMap[tag] < 0, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %" PetscInt_FMT, tag); 11230598e1a2SLisandro Dalcin gmsh->nodeMap[tag] = n; 11240598e1a2SLisandro Dalcin } 11250598e1a2SLisandro Dalcin } 11260598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11270598e1a2SLisandro Dalcin } 11280598e1a2SLisandro Dalcin 11299371c9d4SSatish Balay static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh) { 11300598e1a2SLisandro Dalcin PetscFunctionBegin; 11310598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11329566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadElements_v41(gmsh, mesh)); break; 11339566063dSJacob Faibussowitsch case 40: PetscCall(GmshReadElements_v40(gmsh, mesh)); break; 11349566063dSJacob Faibussowitsch default: PetscCall(GmshReadElements_v22(gmsh, mesh)); break; 11350598e1a2SLisandro Dalcin } 11360598e1a2SLisandro Dalcin 11370598e1a2SLisandro Dalcin { /* Reorder elements by codimension and polytope type */ 11380598e1a2SLisandro Dalcin PetscInt ne = mesh->numElems; 11390598e1a2SLisandro Dalcin GmshElement *elements = mesh->elements; 1140066ea43fSLisandro Dalcin PetscInt keymap[GMSH_NUM_POLYTOPES], nk = 0; 1141066ea43fSLisandro Dalcin PetscInt offset[GMSH_NUM_POLYTOPES + 1], e, k; 11420598e1a2SLisandro Dalcin 1143066ea43fSLisandro Dalcin for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT; 11449566063dSJacob Faibussowitsch PetscCall(PetscMemzero(offset, sizeof(offset))); 11450598e1a2SLisandro Dalcin 1146066ea43fSLisandro Dalcin keymap[GMSH_TET] = nk++; 1147066ea43fSLisandro Dalcin keymap[GMSH_HEX] = nk++; 1148066ea43fSLisandro Dalcin keymap[GMSH_PRI] = nk++; 1149066ea43fSLisandro Dalcin keymap[GMSH_PYR] = nk++; 1150066ea43fSLisandro Dalcin keymap[GMSH_TRI] = nk++; 1151066ea43fSLisandro Dalcin keymap[GMSH_QUA] = nk++; 1152066ea43fSLisandro Dalcin keymap[GMSH_SEG] = nk++; 1153066ea43fSLisandro Dalcin keymap[GMSH_VTX] = nk++; 11540598e1a2SLisandro Dalcin 11559566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(mesh->numElems, &mesh->elements)); 11560598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope] 11570598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) offset[1 + key(e)]++; 11580598e1a2SLisandro Dalcin for (k = 1; k < nk; ++k) offset[k] += offset[k - 1]; 11590598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e]; 11600598e1a2SLisandro Dalcin #undef key 11619566063dSJacob Faibussowitsch PetscCall(GmshElementsDestroy(&elements)); 1162066ea43fSLisandro Dalcin } 11630598e1a2SLisandro Dalcin 1164066ea43fSLisandro Dalcin { /* Mesh dimension and order */ 1165066ea43fSLisandro Dalcin GmshElement *elem = mesh->numElems ? mesh->elements : NULL; 1166066ea43fSLisandro Dalcin mesh->dim = elem ? GmshCellMap[elem->cellType].dim : 0; 1167066ea43fSLisandro Dalcin mesh->order = elem ? GmshCellMap[elem->cellType].order : 0; 11680598e1a2SLisandro Dalcin } 11690598e1a2SLisandro Dalcin 11700598e1a2SLisandro Dalcin { 11710598e1a2SLisandro Dalcin PetscBT vtx; 11720598e1a2SLisandro Dalcin PetscInt dim = mesh->dim, e, n, v; 11730598e1a2SLisandro Dalcin 11749566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(mesh->numNodes, &vtx)); 11750598e1a2SLisandro Dalcin 11760598e1a2SLisandro Dalcin /* Compute number of cells and set of vertices */ 11770598e1a2SLisandro Dalcin mesh->numCells = 0; 11780598e1a2SLisandro Dalcin for (e = 0; e < mesh->numElems; ++e) { 11790598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 11800598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) mesh->numCells++; 1181*48a46eb9SPierre Jolivet for (v = 0; v < elem->numVerts; v++) PetscCall(PetscBTSet(vtx, elem->nodes[v])); 11820598e1a2SLisandro Dalcin } 11830598e1a2SLisandro Dalcin 11840598e1a2SLisandro Dalcin /* Compute numbering for vertices */ 11850598e1a2SLisandro Dalcin mesh->numVerts = 0; 11869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mesh->numNodes, &mesh->vertexMap)); 11879371c9d4SSatish Balay for (n = 0; n < mesh->numNodes; ++n) mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT; 11880598e1a2SLisandro Dalcin 11899566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&vtx)); 11900598e1a2SLisandro Dalcin } 11910598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11920598e1a2SLisandro Dalcin } 11930598e1a2SLisandro Dalcin 11949371c9d4SSatish Balay static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh) { 11950598e1a2SLisandro Dalcin PetscInt n; 11960598e1a2SLisandro Dalcin 11970598e1a2SLisandro Dalcin PetscFunctionBegin; 11989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mesh->numNodes, &mesh->periodMap)); 11990598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n; 12000598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 12019566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadPeriodic_v41(gmsh, mesh->periodMap)); break; 12029566063dSJacob Faibussowitsch default: PetscCall(GmshReadPeriodic_v40(gmsh, mesh->periodMap)); break; 12030598e1a2SLisandro Dalcin } 12040598e1a2SLisandro Dalcin 12059dddd249SSatish Balay /* Find canonical primary nodes */ 12060598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12079371c9d4SSatish Balay while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]]) mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]]; 12080598e1a2SLisandro Dalcin 12099dddd249SSatish Balay /* Renumber vertices (filter out correspondings) */ 12100598e1a2SLisandro Dalcin mesh->numVerts = 0; 12110598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12120598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12139dddd249SSatish Balay if (mesh->periodMap[n] == n) /* is primary */ 12140598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->numVerts++; 12150598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12160598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12179dddd249SSatish Balay if (mesh->periodMap[n] != n) /* is corresponding */ 12180598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]]; 12190598e1a2SLisandro Dalcin PetscFunctionReturn(0); 12200598e1a2SLisandro Dalcin } 12210598e1a2SLisandro Dalcin 1222066ea43fSLisandro Dalcin #define DM_POLYTOPE_VERTEX DM_POLYTOPE_POINT 1223066ea43fSLisandro Dalcin static const DMPolytopeType DMPolytopeMap[] = { 1224066ea43fSLisandro Dalcin /* GMSH_VTX */ DM_POLYTOPE_VERTEX, 1225066ea43fSLisandro Dalcin /* GMSH_SEG */ DM_POLYTOPE_SEGMENT, 1226066ea43fSLisandro Dalcin /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE, 1227066ea43fSLisandro Dalcin /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL, 1228066ea43fSLisandro Dalcin /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON, 1229066ea43fSLisandro Dalcin /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON, 1230066ea43fSLisandro Dalcin /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM, 12319371c9d4SSatish Balay /* GMSH_PYR */ DM_POLYTOPE_PYRAMID, DM_POLYTOPE_UNKNOWN}; 12320598e1a2SLisandro Dalcin 12339371c9d4SSatish Balay static inline DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType) { 1234066ea43fSLisandro Dalcin return DMPolytopeMap[GmshCellMap[cellType].polytope]; 1235066ea43fSLisandro Dalcin } 1236066ea43fSLisandro Dalcin 12379371c9d4SSatish Balay static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem) { 1238066ea43fSLisandro Dalcin DM K; 1239066ea43fSLisandro Dalcin PetscSpace P; 1240066ea43fSLisandro Dalcin PetscDualSpace Q; 1241066ea43fSLisandro Dalcin PetscQuadrature q, fq; 1242066ea43fSLisandro Dalcin PetscBool isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 1243066ea43fSLisandro Dalcin PetscBool endpoint = PETSC_TRUE; 1244066ea43fSLisandro Dalcin char name[32]; 1245066ea43fSLisandro Dalcin 1246066ea43fSLisandro Dalcin PetscFunctionBegin; 1247066ea43fSLisandro Dalcin /* Create space */ 12489566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(comm, &P)); 12499566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 12509566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(P, isTensor)); 12519566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(P, Nc)); 12529566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(P, dim)); 12539566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(P, k, PETSC_DETERMINE)); 1254066ea43fSLisandro Dalcin if (prefix) { 12559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, prefix)); 12569566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetFromOptions(P)); 12579566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)P, NULL)); 12589566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &k, NULL)); 1259066ea43fSLisandro Dalcin } 12609566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(P)); 1261066ea43fSLisandro Dalcin /* Create dual space */ 12629566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(comm, &Q)); 12639566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 12649566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor)); 12659566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity)); 12669566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0)); 12679566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 12689566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(Q, k)); 12699566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 12709566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(Q, K)); 12719566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 1272066ea43fSLisandro Dalcin if (prefix) { 12739566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, prefix)); 12749566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFromOptions(Q)); 12759566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)Q, NULL)); 1276066ea43fSLisandro Dalcin } 12779566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(Q)); 1278066ea43fSLisandro Dalcin /* Create quadrature */ 1279066ea43fSLisandro Dalcin if (isSimplex) { 12809566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, k + 1, -1, +1, &q)); 12819566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim - 1, 1, k + 1, -1, +1, &fq)); 1282066ea43fSLisandro Dalcin } else { 12839566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, k + 1, -1, +1, &q)); 12849566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim - 1, 1, k + 1, -1, +1, &fq)); 1285066ea43fSLisandro Dalcin } 1286066ea43fSLisandro Dalcin /* Create finite element */ 12879566063dSJacob Faibussowitsch PetscCall(PetscFECreate(comm, fem)); 12889566063dSJacob Faibussowitsch PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 12899566063dSJacob Faibussowitsch PetscCall(PetscFESetNumComponents(*fem, Nc)); 12909566063dSJacob Faibussowitsch PetscCall(PetscFESetBasisSpace(*fem, P)); 12919566063dSJacob Faibussowitsch PetscCall(PetscFESetDualSpace(*fem, Q)); 12929566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(*fem, q)); 12939566063dSJacob Faibussowitsch PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 1294066ea43fSLisandro Dalcin if (prefix) { 12959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, prefix)); 12969566063dSJacob Faibussowitsch PetscCall(PetscFESetFromOptions(*fem)); 12979566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*fem, NULL)); 1298066ea43fSLisandro Dalcin } 12999566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(*fem)); 1300066ea43fSLisandro Dalcin /* Cleanup */ 13019566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&P)); 13029566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&Q)); 13039566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 13049566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&fq)); 1305066ea43fSLisandro Dalcin /* Set finite element name */ 130663a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, sizeof(name), "%s%" PetscInt_FMT, isSimplex ? "P" : "Q", k)); 13079566063dSJacob Faibussowitsch PetscCall(PetscFESetName(*fem, name)); 1308066ea43fSLisandro Dalcin PetscFunctionReturn(0); 130996ca5757SLisandro Dalcin } 131096ca5757SLisandro Dalcin 1311d6689ca9SLisandro Dalcin /*@C 1312d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 1313d6689ca9SLisandro Dalcin 1314d6689ca9SLisandro Dalcin + comm - The MPI communicator 1315d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 1316d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 1317d6689ca9SLisandro Dalcin 1318d6689ca9SLisandro Dalcin Output Parameter: 1319d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 1320d6689ca9SLisandro Dalcin 1321d6689ca9SLisandro Dalcin Level: beginner 1322d6689ca9SLisandro Dalcin 1323db781477SPatrick Sanan .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 1324d6689ca9SLisandro Dalcin @*/ 13259371c9d4SSatish Balay PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) { 1326d6689ca9SLisandro Dalcin PetscViewer viewer; 1327d6689ca9SLisandro Dalcin PetscMPIInt rank; 1328d6689ca9SLisandro Dalcin int fileType; 1329d6689ca9SLisandro Dalcin PetscViewerType vtype; 1330d6689ca9SLisandro Dalcin 1331d6689ca9SLisandro Dalcin PetscFunctionBegin; 13329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1333d6689ca9SLisandro Dalcin 1334d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 1335dd400576SPatrick Sanan if (rank == 0) { 13360598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1337d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1338d6689ca9SLisandro Dalcin int snum; 1339d6689ca9SLisandro Dalcin float version; 1340a8d4e440SJunchao Zhang int fileFormat; 1341d6689ca9SLisandro Dalcin 13429566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gmsh, 1)); 13439566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer)); 13449566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII)); 13459566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ)); 13469566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(gmsh->viewer, filename)); 1347d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 13489566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 13499566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$MeshFormat", line)); 13509566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 2)); 1351d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 1352a8d4e440SJunchao Zhang fileFormat = (int)roundf(version * 10); 135308401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 1354a8d4e440SJunchao Zhang PetscCheck(fileFormat >= 22, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 13551dca8a05SBarry Smith PetscCheck((int)version != 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 1356a8d4e440SJunchao Zhang PetscCheck(fileFormat <= 41, PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 13579566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&gmsh->viewer)); 1358d6689ca9SLisandro Dalcin } 13599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&fileType, 1, MPI_INT, 0, comm)); 1360d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 1361d6689ca9SLisandro Dalcin 1362d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 13639566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 13649566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, vtype)); 13659566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 13669566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 13679566063dSJacob Faibussowitsch PetscCall(DMPlexCreateGmsh(comm, viewer, interpolate, dm)); 13689566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1369d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 1370d6689ca9SLisandro Dalcin } 1371d6689ca9SLisandro Dalcin 1372331830f3SMatthew G. Knepley /*@ 13737d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 1374331830f3SMatthew G. Knepley 1375d083f849SBarry Smith Collective 1376331830f3SMatthew G. Knepley 1377331830f3SMatthew G. Knepley Input Parameters: 1378331830f3SMatthew G. Knepley + comm - The MPI communicator 1379331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1380331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1381331830f3SMatthew G. Knepley 1382331830f3SMatthew G. Knepley Output Parameter: 1383331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1384331830f3SMatthew G. Knepley 1385df93260fSMatthew G. Knepley Options Database: 1386df93260fSMatthew G. Knepley + -dm_plex_gmsh_hybrid - Force triangular prisms to use tensor order 1387df93260fSMatthew G. Knepley . -dm_plex_gmsh_periodic - Read Gmsh periodic section and construct a periodic Plex 1388df93260fSMatthew G. Knepley . -dm_plex_gmsh_highorder - Generate high-order coordinates 1389df93260fSMatthew G. Knepley . -dm_plex_gmsh_project - Project high-order coordinates to a different space, use the prefix dm_plex_gmsh_project_ to define the space 1390df93260fSMatthew G. Knepley . -dm_plex_gmsh_use_regions - Generate labels with region names 1391df93260fSMatthew G. Knepley . -dm_plex_gmsh_mark_vertices - Add vertices to generated labels 1392df93260fSMatthew G. Knepley . -dm_plex_gmsh_multiple_tags - Allow multiple tags for default labels 1393df93260fSMatthew G. Knepley - -dm_plex_gmsh_spacedim <d> - Embedding space dimension, if different from topological dimension 1394df93260fSMatthew G. Knepley 1395df93260fSMatthew G. Knepley Note: 1396df93260fSMatthew G. Knepley The Gmsh file format is described in http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1397df93260fSMatthew G. Knepley 1398df93260fSMatthew G. Knepley 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. Instead, -dm_plex_gmsh_use_regions creates labels based on the region names from the PhysicalNames section, and all tags are used. 1399331830f3SMatthew G. Knepley 1400331830f3SMatthew G. Knepley Level: beginner 1401331830f3SMatthew G. Knepley 1402db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()` 1403331830f3SMatthew G. Knepley @*/ 14049371c9d4SSatish Balay PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) { 14050598e1a2SLisandro Dalcin GmshMesh *mesh = NULL; 140611c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 14070598e1a2SLisandro Dalcin PetscBT periodicVerts = NULL; 14080598e1a2SLisandro Dalcin PetscBT periodicCells = NULL; 14096858538eSMatthew G. Knepley DM cdm, cdmCell = NULL; 14106858538eSMatthew G. Knepley PetscSection cs, csCell = NULL; 14116858538eSMatthew G. Knepley Vec coordinates, coordinatesCell; 1412a45dabc8SMatthew G. Knepley DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL, marker = NULL, *regionSets; 1413066ea43fSLisandro Dalcin PetscInt dim = 0, coordDim = -1, order = 0; 14140598e1a2SLisandro Dalcin PetscInt numNodes = 0, numElems = 0, numVerts = 0, numCells = 0; 1415066ea43fSLisandro Dalcin PetscInt cell, cone[8], e, n, v, d; 1416df93260fSMatthew G. Knepley PetscBool binary, useregions = PETSC_FALSE, markvertices = PETSC_FALSE, multipleTags = PETSC_FALSE; 14170598e1a2SLisandro Dalcin PetscBool hybrid = interpolate, periodic = PETSC_TRUE; 1418066ea43fSLisandro Dalcin PetscBool highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE; 1419066ea43fSLisandro Dalcin PetscBool isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE; 14200598e1a2SLisandro Dalcin PetscMPIInt rank; 1421331830f3SMatthew G. Knepley 1422331830f3SMatthew G. Knepley PetscFunctionBegin; 14239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1424d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)viewer); 1425d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMPlex Gmsh options"); 1426df93260fSMatthew G. Knepley PetscCall(PetscOptionsBool("-dm_plex_gmsh_hybrid", "Force triangular prisms to use tensor order", "DMPlexCreateGmsh", hybrid, &hybrid, NULL)); 14279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_periodic", "Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL)); 14289566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_highorder", "Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet)); 14299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL)); 14309566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_regions", "Generate labels with region names", "DMPlexCreateGmsh", useregions, &useregions, NULL)); 14319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_mark_vertices", "Add vertices to generated labels", "DMPlexCreateGmsh", markvertices, &markvertices, NULL)); 1432df93260fSMatthew G. Knepley PetscCall(PetscOptionsBool("-dm_plex_gmsh_multiple_tags", "Allow multiple tags for default labels", "DMPlexCreateGmsh", multipleTags, &multipleTags, NULL)); 14339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE)); 1434d0609cedSBarry Smith PetscOptionsHeadEnd(); 1435d0609cedSBarry Smith PetscOptionsEnd(); 14360598e1a2SLisandro Dalcin 14379566063dSJacob Faibussowitsch PetscCall(GmshCellInfoSetUp()); 143811c56cb0SLisandro Dalcin 14399566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 14409566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 14419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_CreateGmsh, *dm, NULL, NULL, NULL)); 144211c56cb0SLisandro Dalcin 14439566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary)); 144411c56cb0SLisandro Dalcin 144511c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 14463b46f529SLisandro Dalcin if (binary) { 144711c56cb0SLisandro Dalcin parentviewer = viewer; 14489566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer)); 144911c56cb0SLisandro Dalcin } 145011c56cb0SLisandro Dalcin 1451dd400576SPatrick Sanan if (rank == 0) { 14520598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1453698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1454698a79b9SLisandro Dalcin PetscBool match; 1455331830f3SMatthew G. Knepley 14569566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gmsh, 1)); 1457698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1458698a79b9SLisandro Dalcin gmsh->binary = binary; 1459698a79b9SLisandro Dalcin 14609566063dSJacob Faibussowitsch PetscCall(GmshMeshCreate(&mesh)); 14610598e1a2SLisandro Dalcin 1462698a79b9SLisandro Dalcin /* Read mesh format */ 14639566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 14649566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$MeshFormat", line)); 14659566063dSJacob Faibussowitsch PetscCall(GmshReadMeshFormat(gmsh)); 14669566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndMeshFormat", line)); 146711c56cb0SLisandro Dalcin 1468dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 14699566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 14709566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$PhysicalNames", line, &match)); 1471dc0ede3bSMatthew G. Knepley if (match) { 14729566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$PhysicalNames", line)); 14739566063dSJacob Faibussowitsch PetscCall(GmshReadPhysicalNames(gmsh, mesh)); 14749566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndPhysicalNames", line)); 1475698a79b9SLisandro Dalcin /* Initial read for entity section */ 14769566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 1477dc0ede3bSMatthew G. Knepley } 147811c56cb0SLisandro Dalcin 1479de78e4feSLisandro Dalcin /* Read entities */ 1480698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 14819566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Entities", line)); 14829566063dSJacob Faibussowitsch PetscCall(GmshReadEntities(gmsh, mesh)); 14839566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndEntities", line)); 1484698a79b9SLisandro Dalcin /* Initial read for nodes section */ 14859566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 1486de78e4feSLisandro Dalcin } 1487de78e4feSLisandro Dalcin 1488698a79b9SLisandro Dalcin /* Read nodes */ 14899566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Nodes", line)); 14909566063dSJacob Faibussowitsch PetscCall(GmshReadNodes(gmsh, mesh)); 14919566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndNodes", line)); 149211c56cb0SLisandro Dalcin 1493698a79b9SLisandro Dalcin /* Read elements */ 14949566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 14959566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Elements", line)); 14969566063dSJacob Faibussowitsch PetscCall(GmshReadElements(gmsh, mesh)); 14979566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndElements", line)); 1498de78e4feSLisandro Dalcin 14990598e1a2SLisandro Dalcin /* Read periodic section (OPTIONAL) */ 1500698a79b9SLisandro Dalcin if (periodic) { 15019566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 15029566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$Periodic", line, &periodic)); 1503ef12879bSLisandro Dalcin } 1504ef12879bSLisandro Dalcin if (periodic) { 15059566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Periodic", line)); 15069566063dSJacob Faibussowitsch PetscCall(GmshReadPeriodic(gmsh, mesh)); 15079566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndPeriodic", line)); 1508698a79b9SLisandro Dalcin } 1509698a79b9SLisandro Dalcin 15109566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->wbuf)); 15119566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->sbuf)); 15129566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->nbuf)); 15130598e1a2SLisandro Dalcin 15140598e1a2SLisandro Dalcin dim = mesh->dim; 1515066ea43fSLisandro Dalcin order = mesh->order; 15160598e1a2SLisandro Dalcin numNodes = mesh->numNodes; 15170598e1a2SLisandro Dalcin numElems = mesh->numElems; 15180598e1a2SLisandro Dalcin numVerts = mesh->numVerts; 15190598e1a2SLisandro Dalcin numCells = mesh->numCells; 1520066ea43fSLisandro Dalcin 1521066ea43fSLisandro Dalcin { 1522066ea43fSLisandro Dalcin GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL; 1523066ea43fSLisandro Dalcin GmshElement *elemB = elemA ? elemA + mesh->numCells - 1 : NULL; 1524066ea43fSLisandro Dalcin int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1; 1525066ea43fSLisandro Dalcin int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1; 1526066ea43fSLisandro Dalcin isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE; 1527066ea43fSLisandro Dalcin isHybrid = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE; 1528066ea43fSLisandro Dalcin hasTetra = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE; 1529066ea43fSLisandro Dalcin } 1530698a79b9SLisandro Dalcin } 1531698a79b9SLisandro Dalcin 1532*48a46eb9SPierre Jolivet if (parentviewer) PetscCall(PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer)); 1533698a79b9SLisandro Dalcin 1534066ea43fSLisandro Dalcin { 1535066ea43fSLisandro Dalcin int buf[6]; 1536698a79b9SLisandro Dalcin 1537066ea43fSLisandro Dalcin buf[0] = (int)dim; 1538066ea43fSLisandro Dalcin buf[1] = (int)order; 1539066ea43fSLisandro Dalcin buf[2] = periodic; 1540066ea43fSLisandro Dalcin buf[3] = isSimplex; 1541066ea43fSLisandro Dalcin buf[4] = isHybrid; 1542066ea43fSLisandro Dalcin buf[5] = hasTetra; 1543066ea43fSLisandro Dalcin 15449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(buf, 6, MPI_INT, 0, comm)); 1545066ea43fSLisandro Dalcin 1546066ea43fSLisandro Dalcin dim = buf[0]; 1547066ea43fSLisandro Dalcin order = buf[1]; 1548066ea43fSLisandro Dalcin periodic = buf[2] ? PETSC_TRUE : PETSC_FALSE; 1549066ea43fSLisandro Dalcin isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE; 1550066ea43fSLisandro Dalcin isHybrid = buf[4] ? PETSC_TRUE : PETSC_FALSE; 1551066ea43fSLisandro Dalcin hasTetra = buf[5] ? PETSC_TRUE : PETSC_FALSE; 1552dea2a3c7SStefano Zampini } 1553dea2a3c7SStefano Zampini 1554066ea43fSLisandro Dalcin if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE; 155508401ef6SPierre Jolivet PetscCheck(!highOrder || !isHybrid, comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet"); 1556066ea43fSLisandro Dalcin 15570598e1a2SLisandro Dalcin /* We do not want this label automatically computed, instead we fill it here */ 15589566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "celltype")); 155911c56cb0SLisandro Dalcin 1560a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 15619566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, numCells + numVerts)); 15620598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 15630598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 15640598e1a2SLisandro Dalcin DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType); 15650598e1a2SLisandro Dalcin if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR; 15669566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(*dm, cell, elem->numVerts)); 15679566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, cell, ctype)); 1568db397164SMichael Lange } 1569*48a46eb9SPierre Jolivet for (v = numCells; v < numCells + numVerts; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 15709566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 157196ca5757SLisandro Dalcin 1572a4bb7517SMichael Lange /* Add cell-vertex connections */ 15730598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 15740598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 15750598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 15760598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 15770598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 15780598e1a2SLisandro Dalcin cone[v] = numCells + vv; 1579db397164SMichael Lange } 15809566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(*dm, cell, cone)); 15819566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, cell, cone)); 1582a4bb7517SMichael Lange } 158396ca5757SLisandro Dalcin 15849566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 15859566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 15869566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1587331830f3SMatthew G. Knepley if (interpolate) { 15885fd9971aSMatthew G. Knepley DM idm; 1589331830f3SMatthew G. Knepley 15909566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 15919566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1592331830f3SMatthew G. Knepley *dm = idm; 1593331830f3SMatthew G. Knepley } 1594917266a4SLisandro Dalcin 1595dd400576SPatrick Sanan if (rank == 0) { 1596a45dabc8SMatthew G. Knepley const PetscInt Nr = useregions ? mesh->numRegions : 0; 159711c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1598d1a54cefSMatthew G. Knepley 15999566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nr, ®ionSets)); 16009566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd)); 16010598e1a2SLisandro Dalcin for (cell = 0, e = 0; e < numElems; ++e) { 16020598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 16036e1dd89aSLawrence Mitchell 16046e1dd89aSLawrence Mitchell /* Create cell sets */ 16050598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) { 16060598e1a2SLisandro Dalcin if (elem->numTags > 0) { 1607df93260fSMatthew G. Knepley PetscInt Nt = elem->numTags, t, r; 1608a45dabc8SMatthew G. Knepley 1609df93260fSMatthew G. Knepley for (t = 0; t < Nt; ++t) { 1610df93260fSMatthew G. Knepley const PetscInt tag = elem->tags[t]; 1611df93260fSMatthew G. Knepley const PetscBool generic = !Nr && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE; 1612df93260fSMatthew G. Knepley 1613df93260fSMatthew G. Knepley if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", cell, tag)); 1614a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1615df93260fSMatthew G. Knepley if (mesh->regionDims[r] != dim) continue; 16169566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], cell, tag)); 1617a45dabc8SMatthew G. Knepley } 16186e1dd89aSLawrence Mitchell } 1619df93260fSMatthew G. Knepley } 1620917266a4SLisandro Dalcin cell++; 16216e1dd89aSLawrence Mitchell } 16220c070f12SSander Arens 16230598e1a2SLisandro Dalcin /* Create face sets */ 16240598e1a2SLisandro Dalcin if (interpolate && elem->dim == dim - 1) { 16250598e1a2SLisandro Dalcin PetscInt joinSize; 16260598e1a2SLisandro Dalcin const PetscInt *join = NULL; 16275ad74b4fSMatthew G. Knepley PetscInt Nt = elem->numTags, t, r; 1628a45dabc8SMatthew G. Knepley 16290598e1a2SLisandro Dalcin /* Find the relevant facet with vertex joins */ 16300598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 16310598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 16320598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 16330598e1a2SLisandro Dalcin cone[v] = vStart + vv; 16340598e1a2SLisandro Dalcin } 16359566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join)); 163663a3b9bcSJacob Faibussowitsch 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); 16375ad74b4fSMatthew G. Knepley for (t = 0; t < Nt; ++t) { 16385ad74b4fSMatthew G. Knepley const PetscInt tag = elem->tags[t]; 1639df93260fSMatthew G. Knepley const PetscBool generic = !Nr && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE; 16405ad74b4fSMatthew G. Knepley 1641df93260fSMatthew G. Knepley if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], tag)); 1642a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1643df93260fSMatthew G. Knepley if (mesh->regionDims[r] != dim - 1) continue; 16449566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], join[0], tag)); 1645a45dabc8SMatthew G. Knepley } 16465ad74b4fSMatthew G. Knepley } 16479566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join)); 16480598e1a2SLisandro Dalcin } 16490598e1a2SLisandro Dalcin 16500c070f12SSander Arens /* Create vertex sets */ 1651df93260fSMatthew G. Knepley if (elem->dim == 0 && markvertices) { 16520598e1a2SLisandro Dalcin if (elem->numTags > 0) { 16530598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[0]; 16540598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1655a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1656a45dabc8SMatthew G. Knepley PetscInt r; 1657a45dabc8SMatthew G. Knepley 16589566063dSJacob Faibussowitsch if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag)); 165981a1af93SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1660df93260fSMatthew G. Knepley if (mesh->regionDims[r] != 0) continue; 16619566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag)); 166281a1af93SMatthew G. Knepley } 166381a1af93SMatthew G. Knepley } 166481a1af93SMatthew G. Knepley } 166581a1af93SMatthew G. Knepley } 166681a1af93SMatthew G. Knepley if (markvertices) { 166781a1af93SMatthew G. Knepley for (v = 0; v < numNodes; ++v) { 166881a1af93SMatthew G. Knepley const PetscInt vv = mesh->vertexMap[v]; 16697d5b9244SMatthew G. Knepley const PetscInt *tags = &mesh->nodelist->tag[v * GMSH_MAX_TAGS]; 16707d5b9244SMatthew G. Knepley PetscInt r, t; 167181a1af93SMatthew G. Knepley 16727d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) { 16737d5b9244SMatthew G. Knepley const PetscInt tag = tags[t]; 1674df93260fSMatthew G. Knepley const PetscBool generic = !Nr && (!t || multipleTags) ? PETSC_TRUE : PETSC_FALSE; 16757d5b9244SMatthew G. Knepley 16767d5b9244SMatthew G. Knepley if (tag == -1) continue; 1677df93260fSMatthew G. Knepley if (generic) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag)); 1678a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 16799566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag)); 16800598e1a2SLisandro Dalcin } 16810598e1a2SLisandro Dalcin } 16820598e1a2SLisandro Dalcin } 16830598e1a2SLisandro Dalcin } 16849566063dSJacob Faibussowitsch PetscCall(PetscFree(regionSets)); 1685a45dabc8SMatthew G. Knepley } 16860598e1a2SLisandro Dalcin 16877dd454faSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 16889371c9d4SSatish Balay enum { 16899371c9d4SSatish Balay n = 4 16909371c9d4SSatish Balay }; 16917dd454faSLisandro Dalcin PetscBool flag[n]; 16927dd454faSLisandro Dalcin 16937dd454faSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 16947dd454faSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 16957dd454faSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 16967dd454faSLisandro Dalcin flag[3] = marker ? PETSC_TRUE : PETSC_FALSE; 16979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 16989566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 16999566063dSJacob Faibussowitsch if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 17009566063dSJacob Faibussowitsch if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 17019566063dSJacob Faibussowitsch if (flag[3]) PetscCall(DMCreateLabel(*dm, "marker")); 17027dd454faSLisandro Dalcin } 17037dd454faSLisandro Dalcin 17040598e1a2SLisandro Dalcin if (periodic) { 17059566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(numVerts, &periodicVerts)); 17060598e1a2SLisandro Dalcin for (n = 0; n < numNodes; ++n) { 17070598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) { 17080598e1a2SLisandro Dalcin if (PetscUnlikely(mesh->periodMap[n] != n)) { 17090598e1a2SLisandro Dalcin PetscInt m = mesh->periodMap[n]; 17109566063dSJacob Faibussowitsch PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[n])); 17119566063dSJacob Faibussowitsch PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[m])); 17120598e1a2SLisandro Dalcin } 17130598e1a2SLisandro Dalcin } 17140598e1a2SLisandro Dalcin } 17159566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(numCells, &periodicCells)); 17160598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 17170598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 17180598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 17190598e1a2SLisandro Dalcin PetscInt nn = elem->nodes[v]; 17200598e1a2SLisandro Dalcin PetscInt vv = mesh->vertexMap[nn]; 17210598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) { 17229371c9d4SSatish Balay PetscCall(PetscBTSet(periodicCells, cell)); 17239371c9d4SSatish Balay break; 17240c070f12SSander Arens } 17250c070f12SSander Arens } 17260c070f12SSander Arens } 172716ad7e67SMichael Lange } 172816ad7e67SMichael Lange 1729066ea43fSLisandro Dalcin /* Setup coordinate DM */ 17300598e1a2SLisandro Dalcin if (coordDim < 0) coordDim = dim; 17319566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, coordDim)); 17329566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1733066ea43fSLisandro Dalcin if (highOrder) { 1734066ea43fSLisandro Dalcin PetscFE fe; 1735066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1736066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_EQUISPACED; 1737066ea43fSLisandro Dalcin 1738066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1739066ea43fSLisandro Dalcin 17409566063dSJacob Faibussowitsch PetscCall(GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe)); 17419566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view")); 17429566063dSJacob Faibussowitsch PetscCall(DMSetField(cdm, 0, NULL, (PetscObject)fe)); 17439566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 17449566063dSJacob Faibussowitsch PetscCall(DMCreateDS(cdm)); 1745066ea43fSLisandro Dalcin } 17466858538eSMatthew G. Knepley if (periodic) { 17476858538eSMatthew G. Knepley PetscCall(DMClone(cdm, &cdmCell)); 17486858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateDM(*dm, cdmCell)); 17496858538eSMatthew G. Knepley } 1750066ea43fSLisandro Dalcin 1751066ea43fSLisandro Dalcin /* Create coordinates */ 1752066ea43fSLisandro Dalcin if (highOrder) { 1753066ea43fSLisandro Dalcin PetscInt maxDof = GmshNumNodes_HEX(order) * coordDim; 1754066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1755066ea43fSLisandro Dalcin PetscSection section; 1756066ea43fSLisandro Dalcin PetscScalar *cellCoords; 1757066ea43fSLisandro Dalcin 17589566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(cdm, NULL)); 17596858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &cs)); 17606858538eSMatthew G. Knepley PetscCall(PetscSectionClone(cs, §ion)); 17619566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(cdm, 0, section)); /* XXX Implement DMPlexSetClosurePermutationLexicographic() */ 1762066ea43fSLisandro Dalcin 17639566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(cdm, &coordinates)); 17649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxDof, &cellCoords)); 1765066ea43fSLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1766066ea43fSLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1767066ea43fSLisandro Dalcin const int *lexorder = GmshCellMap[elem->cellType].lexorder(); 1768b9bf55e5SMatthew G. Knepley int s = 0; 1769066ea43fSLisandro Dalcin for (n = 0; n < elem->numNodes; ++n) { 1770b9bf55e5SMatthew G. Knepley while (lexorder[n + s] < 0) ++s; 1771b9bf55e5SMatthew G. Knepley const PetscInt node = elem->nodes[lexorder[n + s]]; 1772b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[(n + s) * coordDim + d] = (PetscReal)coords[node * 3 + d]; 1773b9bf55e5SMatthew G. Knepley } 1774b9bf55e5SMatthew G. Knepley if (s) { 1775b9bf55e5SMatthew G. Knepley /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */ 1776b9bf55e5SMatthew G. Knepley PetscReal quaCenterWeights[9] = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25}; 1777b9bf55e5SMatthew G. Knepley /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */ 17789371c9d4SSatish Balay 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}; 17799371c9d4SSatish Balay 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}; 17809371c9d4SSatish Balay 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}; 17819371c9d4SSatish Balay 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}; 17829371c9d4SSatish Balay 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}; 17839371c9d4SSatish Balay 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}; 17849371c9d4SSatish Balay 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}; 1785b9bf55e5SMatthew G. Knepley PetscReal *sdWeights2[9] = {NULL, NULL, NULL, NULL, quaCenterWeights, NULL, NULL, NULL, NULL}; 17869371c9d4SSatish Balay PetscReal *sdWeights3[27] = {NULL, NULL, NULL, NULL, hexBottomWeights, NULL, NULL, NULL, NULL, NULL, hexFrontWeights, NULL, hexLeftWeights, hexCenterWeights, hexRightWeights, NULL, hexBackWeights, NULL, 1787b9bf55e5SMatthew G. Knepley NULL, NULL, NULL, NULL, hexTopWeights, NULL, NULL, NULL, NULL}; 1788b9bf55e5SMatthew G. Knepley PetscReal **sdWeights[4] = {NULL, NULL, sdWeights2, sdWeights3}; 1789b9bf55e5SMatthew G. Knepley 1790b9bf55e5SMatthew G. Knepley /* Missing entries in serendipity cell, only works for 8-node quad and 20-node hex */ 1791b9bf55e5SMatthew G. Knepley for (n = 0; n < elem->numNodes + s; ++n) { 1792b9bf55e5SMatthew G. Knepley if (lexorder[n] >= 0) continue; 1793b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[n * coordDim + d] = 0.0; 1794b9bf55e5SMatthew G. Knepley for (int bn = 0; bn < elem->numNodes + s; ++bn) { 1795b9bf55e5SMatthew G. Knepley if (lexorder[bn] < 0) continue; 1796b9bf55e5SMatthew G. Knepley const PetscReal *weights = sdWeights[coordDim][n]; 1797b9bf55e5SMatthew G. Knepley const PetscInt bnode = elem->nodes[lexorder[bn]]; 1798b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[n * coordDim + d] += weights[bn] * (PetscReal)coords[bnode * 3 + d]; 1799b9bf55e5SMatthew G. Knepley } 1800b9bf55e5SMatthew G. Knepley } 1801066ea43fSLisandro Dalcin } 18029566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES)); 1803066ea43fSLisandro Dalcin } 18049566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 18059566063dSJacob Faibussowitsch PetscCall(PetscFree(cellCoords)); 1806066ea43fSLisandro Dalcin } else { 1807066ea43fSLisandro Dalcin PetscInt *nodeMap; 1808066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1809066ea43fSLisandro Dalcin PetscScalar *pointCoords; 1810066ea43fSLisandro Dalcin 18116858538eSMatthew G. Knepley PetscCall(DMGetCoordinateSection(*dm, &cs)); 18126858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(cs, 1)); 18136858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(cs, 0, coordDim)); 18146858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(cs, numCells, numCells + numVerts)); 18156858538eSMatthew G. Knepley for (v = numCells; v < numCells + numVerts; ++v) { 18166858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(cs, v, coordDim)); 18176858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(cs, v, 0, coordDim)); 1818f45c9589SStefano Zampini } 18196858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(cs)); 18206858538eSMatthew G. Knepley 18216858538eSMatthew G. Knepley /* We need to localize coordinates on cells */ 18220598e1a2SLisandro Dalcin if (periodic) { 18236858538eSMatthew G. Knepley PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)cdmCell), &csCell)); 18246858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(csCell, 1)); 18256858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(csCell, 0, coordDim)); 18266858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(csCell, 0, numCells)); 18270598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 18280598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 18290598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 18300598e1a2SLisandro Dalcin PetscInt dof = elem->numVerts * coordDim; 18316858538eSMatthew G. Knepley 18326858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(csCell, cell, dof)); 18336858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(csCell, cell, 0, dof)); 1834f45c9589SStefano Zampini } 1835f45c9589SStefano Zampini } 18366858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(csCell)); 18376858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateSection(*dm, PETSC_DETERMINE, csCell)); 1838f45c9589SStefano Zampini } 1839066ea43fSLisandro Dalcin 18409566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(cdm, &coordinates)); 18419566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &pointCoords)); 18426858538eSMatthew G. Knepley PetscCall(PetscMalloc1(numVerts, &nodeMap)); 18436858538eSMatthew G. Knepley for (n = 0; n < numNodes; n++) 18446858538eSMatthew G. Knepley if (mesh->vertexMap[n] >= 0) nodeMap[mesh->vertexMap[n]] = n; 18456858538eSMatthew G. Knepley for (v = 0; v < numVerts; ++v) { 18466858538eSMatthew G. Knepley PetscInt off, node = nodeMap[v]; 18476858538eSMatthew G. Knepley 18486858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cs, numCells + v, &off)); 18496858538eSMatthew G. Knepley for (d = 0; d < coordDim; ++d) pointCoords[off + d] = (PetscReal)coords[node * 3 + d]; 18506858538eSMatthew G. Knepley } 18516858538eSMatthew G. Knepley PetscCall(VecRestoreArray(coordinates, &pointCoords)); 18526858538eSMatthew G. Knepley PetscCall(PetscFree(nodeMap)); 18536858538eSMatthew G. Knepley 18540598e1a2SLisandro Dalcin if (periodic) { 18556858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdmCell, &coordinatesCell)); 18566858538eSMatthew G. Knepley PetscCall(VecGetArray(coordinatesCell, &pointCoords)); 18570598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 18580598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 18590598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 18600598e1a2SLisandro Dalcin PetscInt off, node; 18616858538eSMatthew G. Knepley for (v = 0; v < elem->numVerts; ++v) cone[v] = elem->nodes[v]; 18626858538eSMatthew G. Knepley PetscCall(DMPlexReorderCell(cdmCell, cell, cone)); 18636858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(csCell, cell, &off)); 18640598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 18659371c9d4SSatish Balay for (node = cone[v], d = 0; d < coordDim; ++d) pointCoords[off++] = (PetscReal)coords[node * 3 + d]; 1866331830f3SMatthew G. Knepley } 1867331830f3SMatthew G. Knepley } 18686858538eSMatthew G. Knepley PetscCall(VecSetBlockSize(coordinatesCell, coordDim)); 18696858538eSMatthew G. Knepley PetscCall(VecRestoreArray(coordinatesCell, &pointCoords)); 18706858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinatesLocal(*dm, coordinatesCell)); 18716858538eSMatthew G. Knepley PetscCall(VecDestroy(&coordinatesCell)); 1872331830f3SMatthew G. Knepley } 18736858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csCell)); 18746858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdmCell)); 1875066ea43fSLisandro Dalcin } 1876066ea43fSLisandro Dalcin 18779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 18789566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, coordDim)); 18799566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 18809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 188111c56cb0SLisandro Dalcin 18829566063dSJacob Faibussowitsch PetscCall(GmshMeshDestroy(&mesh)); 18839566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&periodicVerts)); 18849566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&periodicCells)); 188511c56cb0SLisandro Dalcin 1886066ea43fSLisandro Dalcin if (highOrder && project) { 1887066ea43fSLisandro Dalcin PetscFE fe; 1888066ea43fSLisandro Dalcin const char prefix[] = "dm_plex_gmsh_project_"; 1889066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1890066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI; 1891066ea43fSLisandro Dalcin 1892066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1893066ea43fSLisandro Dalcin 18949566063dSJacob Faibussowitsch PetscCall(GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe)); 18959566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view")); 18969566063dSJacob Faibussowitsch PetscCall(DMProjectCoordinates(*dm, fe)); 18979566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1898066ea43fSLisandro Dalcin } 1899066ea43fSLisandro Dalcin 19009566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_CreateGmsh, *dm, NULL, NULL, NULL)); 1901331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1902331830f3SMatthew G. Knepley } 1903