1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3730356f6SLisandro Dalcin #include <petsc/private/hashmapi.h> 4331830f3SMatthew G. Knepley 5066ea43fSLisandro Dalcin #include <../src/dm/impls/plex/gmshlex.h> 6066ea43fSLisandro Dalcin 7066ea43fSLisandro Dalcin #define GMSH_LEXORDER_ITEM(T, p) \ 8066ea43fSLisandro Dalcin static int *GmshLexOrder_##T##_##p(void) \ 9066ea43fSLisandro Dalcin { \ 10066ea43fSLisandro Dalcin static int Gmsh_LexOrder_##T##_##p[GmshNumNodes_##T(p)] = {-1}; \ 11066ea43fSLisandro Dalcin int *lex = Gmsh_LexOrder_##T##_##p; \ 12066ea43fSLisandro Dalcin if (lex[0] == -1) (void)GmshLexOrder_##T(p, lex, 0); \ 13066ea43fSLisandro Dalcin return lex; \ 14066ea43fSLisandro Dalcin } 15066ea43fSLisandro Dalcin 16b9bf55e5SMatthew G. Knepley static int *GmshLexOrder_QUA_2_Serendipity(void) 17b9bf55e5SMatthew G. Knepley { 18b9bf55e5SMatthew G. Knepley static int Gmsh_LexOrder_QUA_2_Serendipity[9] = {-1}; 19b9bf55e5SMatthew G. Knepley int *lex = Gmsh_LexOrder_QUA_2_Serendipity; 20b9bf55e5SMatthew G. Knepley if (lex[0] == -1) { 21b9bf55e5SMatthew G. Knepley /* Vertices */ 22b9bf55e5SMatthew G. Knepley lex[0] = 0; lex[2] = 1; lex[8] = 2; lex[6] = 3; 23b9bf55e5SMatthew G. Knepley /* Edges */ 24b9bf55e5SMatthew G. Knepley lex[1] = 4; lex[5] = 5; lex[7] = 6; lex[3] = 7; 25b9bf55e5SMatthew G. Knepley /* Cell */ 26b9bf55e5SMatthew G. Knepley lex[4] = -8; 27b9bf55e5SMatthew G. Knepley } 28b9bf55e5SMatthew G. Knepley return lex; 29b9bf55e5SMatthew G. Knepley } 30b9bf55e5SMatthew G. Knepley 31b9bf55e5SMatthew G. Knepley static int *GmshLexOrder_HEX_2_Serendipity(void) 32b9bf55e5SMatthew G. Knepley { 33b9bf55e5SMatthew G. Knepley static int Gmsh_LexOrder_HEX_2_Serendipity[27] = {-1}; 34b9bf55e5SMatthew G. Knepley int *lex = Gmsh_LexOrder_HEX_2_Serendipity; 35b9bf55e5SMatthew G. Knepley if (lex[0] == -1) { 36b9bf55e5SMatthew G. Knepley /* Vertices */ 37b9bf55e5SMatthew G. Knepley lex[ 0] = 0; lex[ 2] = 1; lex[ 8] = 2; lex[ 6] = 3; 38b9bf55e5SMatthew G. Knepley lex[18] = 4; lex[20] = 5; lex[26] = 6; lex[24] = 7; 39b9bf55e5SMatthew G. Knepley /* Edges */ 40b9bf55e5SMatthew G. Knepley lex[ 1] = 8; lex[ 3] = 9; lex[ 9] = 10; lex[ 5] = 11; 41b9bf55e5SMatthew G. Knepley lex[11] = 12; lex[ 7] = 13; lex[17] = 14; lex[15] = 15; 42b9bf55e5SMatthew G. Knepley lex[19] = 16; lex[21] = 17; lex[23] = 18; lex[25] = 19; 43b9bf55e5SMatthew G. Knepley /* Faces */ 44b9bf55e5SMatthew G. Knepley lex[ 4] = -20; lex[10] = -21; lex[12] = -22; lex[14] = -23; 45b9bf55e5SMatthew G. Knepley lex[16] = -24; lex[22] = -25; 46b9bf55e5SMatthew G. Knepley /* Cell */ 47b9bf55e5SMatthew G. Knepley lex[13] = -26; 48b9bf55e5SMatthew G. Knepley } 49b9bf55e5SMatthew G. Knepley return lex; 50b9bf55e5SMatthew G. Knepley } 51b9bf55e5SMatthew G. Knepley 52066ea43fSLisandro Dalcin #define GMSH_LEXORDER_LIST(T) \ 53066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 1) \ 54066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 2) \ 55066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 3) \ 56066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 4) \ 57066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 5) \ 58066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 6) \ 59066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 7) \ 60066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 8) \ 61066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 9) \ 62066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 10) 63066ea43fSLisandro Dalcin 64066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(VTX, 0) 65066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(SEG) 66066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TRI) 67066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(QUA) 68066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TET) 69066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(HEX) 70066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PRI) 71066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PYR) 72066ea43fSLisandro Dalcin 73066ea43fSLisandro Dalcin typedef enum { 74066ea43fSLisandro Dalcin GMSH_VTX = 0, 75066ea43fSLisandro Dalcin GMSH_SEG = 1, 76066ea43fSLisandro Dalcin GMSH_TRI = 2, 77066ea43fSLisandro Dalcin GMSH_QUA = 3, 78066ea43fSLisandro Dalcin GMSH_TET = 4, 79066ea43fSLisandro Dalcin GMSH_HEX = 5, 80066ea43fSLisandro Dalcin GMSH_PRI = 6, 81066ea43fSLisandro Dalcin GMSH_PYR = 7, 82066ea43fSLisandro Dalcin GMSH_NUM_POLYTOPES = 8 83066ea43fSLisandro Dalcin } GmshPolytopeType; 84066ea43fSLisandro Dalcin 85de78e4feSLisandro Dalcin typedef struct { 860598e1a2SLisandro Dalcin int cellType; 87066ea43fSLisandro Dalcin int polytope; 880598e1a2SLisandro Dalcin int dim; 890598e1a2SLisandro Dalcin int order; 90066ea43fSLisandro Dalcin int numVerts; 910598e1a2SLisandro Dalcin int numNodes; 92066ea43fSLisandro Dalcin int* (*lexorder)(void); 930598e1a2SLisandro Dalcin } GmshCellInfo; 940598e1a2SLisandro Dalcin 95066ea43fSLisandro Dalcin #define GmshCellEntry(cellType, polytope, dim, order) \ 96066ea43fSLisandro Dalcin {cellType, GMSH_##polytope, dim, order, \ 97066ea43fSLisandro Dalcin GmshNumNodes_##polytope(1), \ 98066ea43fSLisandro Dalcin GmshNumNodes_##polytope(order), \ 99066ea43fSLisandro Dalcin GmshLexOrder_##polytope##_##order} 1000598e1a2SLisandro Dalcin 1010598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = { 102066ea43fSLisandro Dalcin GmshCellEntry( 15, VTX, 0, 0), 1030598e1a2SLisandro Dalcin 104066ea43fSLisandro Dalcin GmshCellEntry( 1, SEG, 1, 1), 105066ea43fSLisandro Dalcin GmshCellEntry( 8, SEG, 1, 2), 106066ea43fSLisandro Dalcin GmshCellEntry( 26, SEG, 1, 3), 107066ea43fSLisandro Dalcin GmshCellEntry( 27, SEG, 1, 4), 108066ea43fSLisandro Dalcin GmshCellEntry( 28, SEG, 1, 5), 109066ea43fSLisandro Dalcin GmshCellEntry( 62, SEG, 1, 6), 110066ea43fSLisandro Dalcin GmshCellEntry( 63, SEG, 1, 7), 111066ea43fSLisandro Dalcin GmshCellEntry( 64, SEG, 1, 8), 112066ea43fSLisandro Dalcin GmshCellEntry( 65, SEG, 1, 9), 113066ea43fSLisandro Dalcin GmshCellEntry( 66, SEG, 1, 10), 1140598e1a2SLisandro Dalcin 115066ea43fSLisandro Dalcin GmshCellEntry( 2, TRI, 2, 1), 116066ea43fSLisandro Dalcin GmshCellEntry( 9, TRI, 2, 2), 117066ea43fSLisandro Dalcin GmshCellEntry( 21, TRI, 2, 3), 118066ea43fSLisandro Dalcin GmshCellEntry( 23, TRI, 2, 4), 119066ea43fSLisandro Dalcin GmshCellEntry( 25, TRI, 2, 5), 120066ea43fSLisandro Dalcin GmshCellEntry( 42, TRI, 2, 6), 121066ea43fSLisandro Dalcin GmshCellEntry( 43, TRI, 2, 7), 122066ea43fSLisandro Dalcin GmshCellEntry( 44, TRI, 2, 8), 123066ea43fSLisandro Dalcin GmshCellEntry( 45, TRI, 2, 9), 124066ea43fSLisandro Dalcin GmshCellEntry( 46, TRI, 2, 10), 1250598e1a2SLisandro Dalcin 126066ea43fSLisandro Dalcin GmshCellEntry( 3, QUA, 2, 1), 127066ea43fSLisandro Dalcin GmshCellEntry( 10, QUA, 2, 2), 128b9bf55e5SMatthew G. Knepley {16, GMSH_QUA, 2, 2, 4, 8, GmshLexOrder_QUA_2_Serendipity}, 129066ea43fSLisandro Dalcin GmshCellEntry( 36, QUA, 2, 3), 130066ea43fSLisandro Dalcin GmshCellEntry( 37, QUA, 2, 4), 131066ea43fSLisandro Dalcin GmshCellEntry( 38, QUA, 2, 5), 132066ea43fSLisandro Dalcin GmshCellEntry( 47, QUA, 2, 6), 133066ea43fSLisandro Dalcin GmshCellEntry( 48, QUA, 2, 7), 134066ea43fSLisandro Dalcin GmshCellEntry( 49, QUA, 2, 8), 135066ea43fSLisandro Dalcin GmshCellEntry( 50, QUA, 2, 9), 136066ea43fSLisandro Dalcin GmshCellEntry( 51, QUA, 2, 10), 1370598e1a2SLisandro Dalcin 138066ea43fSLisandro Dalcin GmshCellEntry( 4, TET, 3, 1), 139066ea43fSLisandro Dalcin GmshCellEntry( 11, TET, 3, 2), 140066ea43fSLisandro Dalcin GmshCellEntry( 29, TET, 3, 3), 141066ea43fSLisandro Dalcin GmshCellEntry( 30, TET, 3, 4), 142066ea43fSLisandro Dalcin GmshCellEntry( 31, TET, 3, 5), 143066ea43fSLisandro Dalcin GmshCellEntry( 71, TET, 3, 6), 144066ea43fSLisandro Dalcin GmshCellEntry( 72, TET, 3, 7), 145066ea43fSLisandro Dalcin GmshCellEntry( 73, TET, 3, 8), 146066ea43fSLisandro Dalcin GmshCellEntry( 74, TET, 3, 9), 147066ea43fSLisandro Dalcin GmshCellEntry( 75, TET, 3, 10), 1480598e1a2SLisandro Dalcin 149066ea43fSLisandro Dalcin GmshCellEntry( 5, HEX, 3, 1), 150066ea43fSLisandro Dalcin GmshCellEntry( 12, HEX, 3, 2), 151b9bf55e5SMatthew G. Knepley {17, GMSH_HEX, 3, 2, 8, 20, GmshLexOrder_HEX_2_Serendipity}, 152066ea43fSLisandro Dalcin GmshCellEntry( 92, HEX, 3, 3), 153066ea43fSLisandro Dalcin GmshCellEntry( 93, HEX, 3, 4), 154066ea43fSLisandro Dalcin GmshCellEntry( 94, HEX, 3, 5), 155066ea43fSLisandro Dalcin GmshCellEntry( 95, HEX, 3, 6), 156066ea43fSLisandro Dalcin GmshCellEntry( 96, HEX, 3, 7), 157066ea43fSLisandro Dalcin GmshCellEntry( 97, HEX, 3, 8), 158066ea43fSLisandro Dalcin GmshCellEntry( 98, HEX, 3, 9), 159066ea43fSLisandro Dalcin GmshCellEntry( -1, HEX, 3, 10), 1600598e1a2SLisandro Dalcin 161066ea43fSLisandro Dalcin GmshCellEntry( 6, PRI, 3, 1), 162066ea43fSLisandro Dalcin GmshCellEntry( 13, PRI, 3, 2), 163066ea43fSLisandro Dalcin GmshCellEntry( 90, PRI, 3, 3), 164066ea43fSLisandro Dalcin GmshCellEntry( 91, PRI, 3, 4), 165066ea43fSLisandro Dalcin GmshCellEntry(106, PRI, 3, 5), 166066ea43fSLisandro Dalcin GmshCellEntry(107, PRI, 3, 6), 167066ea43fSLisandro Dalcin GmshCellEntry(108, PRI, 3, 7), 168066ea43fSLisandro Dalcin GmshCellEntry(109, PRI, 3, 8), 169066ea43fSLisandro Dalcin GmshCellEntry(110, PRI, 3, 9), 170066ea43fSLisandro Dalcin GmshCellEntry( -1, PRI, 3, 10), 1710598e1a2SLisandro Dalcin 172066ea43fSLisandro Dalcin GmshCellEntry( 7, PYR, 3, 1), 173066ea43fSLisandro Dalcin GmshCellEntry( 14, PYR, 3, 2), 174066ea43fSLisandro Dalcin GmshCellEntry(118, PYR, 3, 3), 175066ea43fSLisandro Dalcin GmshCellEntry(119, PYR, 3, 4), 176066ea43fSLisandro Dalcin GmshCellEntry(120, PYR, 3, 5), 177066ea43fSLisandro Dalcin GmshCellEntry(121, PYR, 3, 6), 178066ea43fSLisandro Dalcin GmshCellEntry(122, PYR, 3, 7), 179066ea43fSLisandro Dalcin GmshCellEntry(123, PYR, 3, 8), 180066ea43fSLisandro Dalcin GmshCellEntry(124, PYR, 3, 9), 181066ea43fSLisandro Dalcin GmshCellEntry( -1, PYR, 3, 10) 1820598e1a2SLisandro Dalcin 1830598e1a2SLisandro Dalcin #if 0 184066ea43fSLisandro Dalcin {20, GMSH_TRI, 2, 3, 3, 9, NULL}, 185066ea43fSLisandro Dalcin {18, GMSH_PRI, 3, 2, 6, 15, NULL}, 186066ea43fSLisandro Dalcin {19, GMSH_PYR, 3, 2, 5, 13, NULL}, 1870598e1a2SLisandro Dalcin #endif 1880598e1a2SLisandro Dalcin }; 1890598e1a2SLisandro Dalcin 1900598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150]; 1910598e1a2SLisandro Dalcin 1920598e1a2SLisandro Dalcin static PetscErrorCode GmshCellInfoSetUp(void) 1930598e1a2SLisandro Dalcin { 1940598e1a2SLisandro Dalcin size_t i, n; 1950598e1a2SLisandro Dalcin static PetscBool called = PETSC_FALSE; 1960598e1a2SLisandro Dalcin 1970598e1a2SLisandro Dalcin if (called) return 0; 1980598e1a2SLisandro Dalcin PetscFunctionBegin; 1990598e1a2SLisandro Dalcin called = PETSC_TRUE; 200dd39110bSPierre Jolivet n = PETSC_STATIC_ARRAY_LENGTH(GmshCellMap); 2010598e1a2SLisandro Dalcin for (i = 0; i < n; ++i) { 2020598e1a2SLisandro Dalcin GmshCellMap[i].cellType = -1; 203066ea43fSLisandro Dalcin GmshCellMap[i].polytope = -1; 2040598e1a2SLisandro Dalcin } 205dd39110bSPierre Jolivet n = PETSC_STATIC_ARRAY_LENGTH(GmshCellTable); 206066ea43fSLisandro Dalcin for (i = 0; i < n; ++i) { 207066ea43fSLisandro Dalcin if (GmshCellTable[i].cellType <= 0) continue; 208066ea43fSLisandro Dalcin GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i]; 209066ea43fSLisandro Dalcin } 2100598e1a2SLisandro Dalcin PetscFunctionReturn(0); 2110598e1a2SLisandro Dalcin } 2120598e1a2SLisandro Dalcin 21394bad497SJacob Faibussowitsch #define GmshCellTypeCheck(ct) PetscMacroReturnStandard( \ 2140598e1a2SLisandro Dalcin const int _ct_ = (int)ct; \ 215dd39110bSPierre Jolivet PetscCheck(_ct_ >= 0 && _ct_ < (int)PETSC_STATIC_ARRAY_LENGTH(GmshCellMap), PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); \ 21694bad497SJacob Faibussowitsch PetscCheck(GmshCellMap[_ct_].cellType == _ct_, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 21794bad497SJacob Faibussowitsch PetscCheck(GmshCellMap[_ct_].polytope != -1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 21894bad497SJacob Faibussowitsch ) 2190598e1a2SLisandro Dalcin 2200598e1a2SLisandro Dalcin typedef struct { 221698a79b9SLisandro Dalcin PetscViewer viewer; 222698a79b9SLisandro Dalcin int fileFormat; 223698a79b9SLisandro Dalcin int dataSize; 224698a79b9SLisandro Dalcin PetscBool binary; 225698a79b9SLisandro Dalcin PetscBool byteSwap; 226698a79b9SLisandro Dalcin size_t wlen; 227698a79b9SLisandro Dalcin void *wbuf; 228698a79b9SLisandro Dalcin size_t slen; 229698a79b9SLisandro Dalcin void *sbuf; 2300598e1a2SLisandro Dalcin PetscInt *nbuf; 2310598e1a2SLisandro Dalcin PetscInt nodeStart; 2320598e1a2SLisandro Dalcin PetscInt nodeEnd; 23333a088b6SMatthew G. Knepley PetscInt *nodeMap; 234698a79b9SLisandro Dalcin } GmshFile; 235de78e4feSLisandro Dalcin 236698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) 237de78e4feSLisandro Dalcin { 238de78e4feSLisandro Dalcin size_t size = count * eltsize; 23911c56cb0SLisandro Dalcin 24011c56cb0SLisandro Dalcin PetscFunctionBegin; 241698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->wbuf)); 2439566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size, &gmsh->wbuf)); 244698a79b9SLisandro Dalcin gmsh->wlen = size; 245de78e4feSLisandro Dalcin } 246698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->wbuf : NULL; 247de78e4feSLisandro Dalcin PetscFunctionReturn(0); 248de78e4feSLisandro Dalcin } 249de78e4feSLisandro Dalcin 250698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) 251698a79b9SLisandro Dalcin { 252698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 253698a79b9SLisandro Dalcin size_t size = count * dataSize; 254698a79b9SLisandro Dalcin 255698a79b9SLisandro Dalcin PetscFunctionBegin; 256698a79b9SLisandro Dalcin if (gmsh->slen < size) { 2579566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->sbuf)); 2589566063dSJacob Faibussowitsch PetscCall(PetscMalloc(size, &gmsh->sbuf)); 259698a79b9SLisandro Dalcin gmsh->slen = size; 260698a79b9SLisandro Dalcin } 261698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->sbuf : NULL; 262698a79b9SLisandro Dalcin PetscFunctionReturn(0); 263698a79b9SLisandro Dalcin } 264698a79b9SLisandro Dalcin 265698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) 266de78e4feSLisandro Dalcin { 267de78e4feSLisandro Dalcin PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype)); 2699566063dSJacob Faibussowitsch if (gmsh->byteSwap) PetscCall(PetscByteSwap(buf, dtype, count)); 270698a79b9SLisandro Dalcin PetscFunctionReturn(0); 271698a79b9SLisandro Dalcin } 272698a79b9SLisandro Dalcin 273698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) 274698a79b9SLisandro Dalcin { 275698a79b9SLisandro Dalcin PetscFunctionBegin; 2769566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING)); 277698a79b9SLisandro Dalcin PetscFunctionReturn(0); 278698a79b9SLisandro Dalcin } 279698a79b9SLisandro Dalcin 280d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) 281d6689ca9SLisandro Dalcin { 282d6689ca9SLisandro Dalcin PetscFunctionBegin; 2839566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(line, Section, match)); 284d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 285d6689ca9SLisandro Dalcin } 286d6689ca9SLisandro Dalcin 287d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) 288d6689ca9SLisandro Dalcin { 289d6689ca9SLisandro Dalcin PetscBool match; 290d6689ca9SLisandro Dalcin 291d6689ca9SLisandro Dalcin PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, Section, line, &match)); 29328b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s",Section); 294d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 295d6689ca9SLisandro Dalcin } 296d6689ca9SLisandro Dalcin 297d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) 298d6689ca9SLisandro Dalcin { 299d6689ca9SLisandro Dalcin PetscBool match; 300d6689ca9SLisandro Dalcin 301d6689ca9SLisandro Dalcin PetscFunctionBegin; 302d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 3039566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 3049566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$Comments", line, &match)); 305d6689ca9SLisandro Dalcin if (!match) break; 306d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 3079566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 3089566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$EndComments", line, &match)); 309d6689ca9SLisandro Dalcin if (match) break; 310d6689ca9SLisandro Dalcin } 311d6689ca9SLisandro Dalcin } 312d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 313d6689ca9SLisandro Dalcin } 314d6689ca9SLisandro Dalcin 315d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) 316d6689ca9SLisandro Dalcin { 317d6689ca9SLisandro Dalcin PetscFunctionBegin; 3189566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 3199566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, EndSection, line)); 320d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 321d6689ca9SLisandro Dalcin } 322d6689ca9SLisandro Dalcin 323698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) 324698a79b9SLisandro Dalcin { 325698a79b9SLisandro Dalcin PetscInt i; 326698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 327698a79b9SLisandro Dalcin 328698a79b9SLisandro Dalcin PetscFunctionBegin; 329698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 3309566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_INT)); 331698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 332698a79b9SLisandro Dalcin int *ibuf = NULL; 3339566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 3349566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_ENUM)); 335698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 336698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 337698a79b9SLisandro Dalcin long *ibuf = NULL; 3389566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 3399566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_LONG)); 340698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 341698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 342698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 3439566063dSJacob Faibussowitsch PetscCall(GmshBufferSizeGet(gmsh, count, &ibuf)); 3449566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, ibuf, count, PETSC_INT64)); 345698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 346698a79b9SLisandro Dalcin } 347698a79b9SLisandro Dalcin PetscFunctionReturn(0); 348698a79b9SLisandro Dalcin } 349698a79b9SLisandro Dalcin 350698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) 351698a79b9SLisandro Dalcin { 352698a79b9SLisandro Dalcin PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_ENUM)); 354698a79b9SLisandro Dalcin PetscFunctionReturn(0); 355698a79b9SLisandro Dalcin } 356698a79b9SLisandro Dalcin 357698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) 358698a79b9SLisandro Dalcin { 359698a79b9SLisandro Dalcin PetscFunctionBegin; 3609566063dSJacob Faibussowitsch PetscCall(GmshRead(gmsh, buf, count, PETSC_DOUBLE)); 361de78e4feSLisandro Dalcin PetscFunctionReturn(0); 362de78e4feSLisandro Dalcin } 363de78e4feSLisandro Dalcin 3649c0edc59SMatthew G. Knepley #define GMSH_MAX_TAGS 16 3657d5b9244SMatthew G. Knepley 366de78e4feSLisandro Dalcin typedef struct { 3670598e1a2SLisandro Dalcin PetscInt id; /* Entity ID */ 3680598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 369de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 370de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 3717d5b9244SMatthew G. Knepley int tags[GMSH_MAX_TAGS]; /* Tag array */ 372de78e4feSLisandro Dalcin } GmshEntity; 373de78e4feSLisandro Dalcin 374de78e4feSLisandro Dalcin typedef struct { 375730356f6SLisandro Dalcin GmshEntity *entity[4]; 376730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 377730356f6SLisandro Dalcin } GmshEntities; 378730356f6SLisandro Dalcin 379698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) 380730356f6SLisandro Dalcin { 381730356f6SLisandro Dalcin PetscInt dim; 382730356f6SLisandro Dalcin 383730356f6SLisandro Dalcin PetscFunctionBegin; 3849566063dSJacob Faibussowitsch PetscCall(PetscNew(entities)); 385730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 3869566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(count[dim], &(*entities)->entity[dim])); 3879566063dSJacob Faibussowitsch PetscCall(PetscHMapICreate(&(*entities)->entityMap[dim])); 388730356f6SLisandro Dalcin } 389730356f6SLisandro Dalcin PetscFunctionReturn(0); 390730356f6SLisandro Dalcin } 391730356f6SLisandro Dalcin 3920598e1a2SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 3930598e1a2SLisandro Dalcin { 3940598e1a2SLisandro Dalcin PetscInt dim; 3950598e1a2SLisandro Dalcin 3960598e1a2SLisandro Dalcin PetscFunctionBegin; 3970598e1a2SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 3980598e1a2SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 3999566063dSJacob Faibussowitsch PetscCall(PetscFree((*entities)->entity[dim])); 4009566063dSJacob Faibussowitsch PetscCall(PetscHMapIDestroy(&(*entities)->entityMap[dim])); 4010598e1a2SLisandro Dalcin } 4029566063dSJacob Faibussowitsch PetscCall(PetscFree((*entities))); 4030598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4040598e1a2SLisandro Dalcin } 4050598e1a2SLisandro Dalcin 406730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 407730356f6SLisandro Dalcin { 408730356f6SLisandro Dalcin PetscFunctionBegin; 4099566063dSJacob Faibussowitsch PetscCall(PetscHMapISet(entities->entityMap[dim], eid, index)); 410730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 411730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 412730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 413730356f6SLisandro Dalcin PetscFunctionReturn(0); 414730356f6SLisandro Dalcin } 415730356f6SLisandro Dalcin 416730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 417730356f6SLisandro Dalcin { 418730356f6SLisandro Dalcin PetscInt index; 419730356f6SLisandro Dalcin 420730356f6SLisandro Dalcin PetscFunctionBegin; 4219566063dSJacob Faibussowitsch PetscCall(PetscHMapIGet(entities->entityMap[dim], eid, &index)); 422730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 423730356f6SLisandro Dalcin PetscFunctionReturn(0); 424730356f6SLisandro Dalcin } 425730356f6SLisandro Dalcin 4260598e1a2SLisandro Dalcin typedef struct { 4270598e1a2SLisandro Dalcin PetscInt *id; /* Node IDs */ 4280598e1a2SLisandro Dalcin double *xyz; /* Coordinates */ 42981a1af93SMatthew G. Knepley PetscInt *tag; /* Physical tag */ 4300598e1a2SLisandro Dalcin } GmshNodes; 4310598e1a2SLisandro Dalcin 4320598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes) 433730356f6SLisandro Dalcin { 434730356f6SLisandro Dalcin PetscFunctionBegin; 4359566063dSJacob Faibussowitsch PetscCall(PetscNew(nodes)); 4369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count*1, &(*nodes)->id)); 4379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count*3, &(*nodes)->xyz)); 4387d5b9244SMatthew G. Knepley PetscCall(PetscMalloc1(count*GMSH_MAX_TAGS, &(*nodes)->tag)); 4390598e1a2SLisandro Dalcin PetscFunctionReturn(0); 440730356f6SLisandro Dalcin } 4410598e1a2SLisandro Dalcin 4420598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes) 4430598e1a2SLisandro Dalcin { 4440598e1a2SLisandro Dalcin PetscFunctionBegin; 4450598e1a2SLisandro Dalcin if (!*nodes) PetscFunctionReturn(0); 4469566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->id)); 4479566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->xyz)); 4489566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes)->tag)); 4499566063dSJacob Faibussowitsch PetscCall(PetscFree((*nodes))); 450730356f6SLisandro Dalcin PetscFunctionReturn(0); 451730356f6SLisandro Dalcin } 452730356f6SLisandro Dalcin 453730356f6SLisandro Dalcin typedef struct { 4540598e1a2SLisandro Dalcin PetscInt id; /* Element ID */ 4550598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 456de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 4570598e1a2SLisandro Dalcin PetscInt numVerts; /* Size of vertex array */ 458de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 4590598e1a2SLisandro Dalcin PetscInt *nodes; /* Vertex/Node array */ 46081a1af93SMatthew G. Knepley PetscInt numTags; /* Size of physical tag array */ 4617d5b9244SMatthew G. Knepley int tags[GMSH_MAX_TAGS]; /* Physical tag array */ 462de78e4feSLisandro Dalcin } GmshElement; 463de78e4feSLisandro Dalcin 4640598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements) 4650598e1a2SLisandro Dalcin { 4660598e1a2SLisandro Dalcin PetscFunctionBegin; 4679566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(count, elements)); 4680598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4690598e1a2SLisandro Dalcin } 4700598e1a2SLisandro Dalcin 4710598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsDestroy(GmshElement **elements) 4720598e1a2SLisandro Dalcin { 4730598e1a2SLisandro Dalcin PetscFunctionBegin; 4740598e1a2SLisandro Dalcin if (!*elements) PetscFunctionReturn(0); 4759566063dSJacob Faibussowitsch PetscCall(PetscFree(*elements)); 4760598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4770598e1a2SLisandro Dalcin } 4780598e1a2SLisandro Dalcin 4790598e1a2SLisandro Dalcin typedef struct { 4800598e1a2SLisandro Dalcin PetscInt dim; 481066ea43fSLisandro Dalcin PetscInt order; 4820598e1a2SLisandro Dalcin GmshEntities *entities; 4830598e1a2SLisandro Dalcin PetscInt numNodes; 4840598e1a2SLisandro Dalcin GmshNodes *nodelist; 4850598e1a2SLisandro Dalcin PetscInt numElems; 4860598e1a2SLisandro Dalcin GmshElement *elements; 4870598e1a2SLisandro Dalcin PetscInt numVerts; 4880598e1a2SLisandro Dalcin PetscInt numCells; 4890598e1a2SLisandro Dalcin PetscInt *periodMap; 4900598e1a2SLisandro Dalcin PetscInt *vertexMap; 4910598e1a2SLisandro Dalcin PetscSegBuffer segbuf; 492a45dabc8SMatthew G. Knepley PetscInt numRegions; 493*90d778b8SLisandro Dalcin PetscInt *regionDims; 494a45dabc8SMatthew G. Knepley PetscInt *regionTags; 495a45dabc8SMatthew G. Knepley char **regionNames; 4960598e1a2SLisandro Dalcin } GmshMesh; 4970598e1a2SLisandro Dalcin 4980598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshCreate(GmshMesh **mesh) 4990598e1a2SLisandro Dalcin { 5000598e1a2SLisandro Dalcin PetscFunctionBegin; 5019566063dSJacob Faibussowitsch PetscCall(PetscNew(mesh)); 5029566063dSJacob Faibussowitsch PetscCall(PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf)); 5030598e1a2SLisandro Dalcin PetscFunctionReturn(0); 5040598e1a2SLisandro Dalcin } 5050598e1a2SLisandro Dalcin 5060598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh) 5070598e1a2SLisandro Dalcin { 508a45dabc8SMatthew G. Knepley PetscInt r; 5090598e1a2SLisandro Dalcin 5100598e1a2SLisandro Dalcin PetscFunctionBegin; 5110598e1a2SLisandro Dalcin if (!*mesh) PetscFunctionReturn(0); 5129566063dSJacob Faibussowitsch PetscCall(GmshEntitiesDestroy(&(*mesh)->entities)); 5139566063dSJacob Faibussowitsch PetscCall(GmshNodesDestroy(&(*mesh)->nodelist)); 5149566063dSJacob Faibussowitsch PetscCall(GmshElementsDestroy(&(*mesh)->elements)); 5159566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh)->periodMap)); 5169566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh)->vertexMap)); 5179566063dSJacob Faibussowitsch PetscCall(PetscSegBufferDestroy(&(*mesh)->segbuf)); 5189566063dSJacob Faibussowitsch for (r = 0; r < (*mesh)->numRegions; ++r) PetscCall(PetscFree((*mesh)->regionNames[r])); 519*90d778b8SLisandro Dalcin PetscCall(PetscFree3((*mesh)->regionDims, (*mesh)->regionTags, (*mesh)->regionNames)); 5209566063dSJacob Faibussowitsch PetscCall(PetscFree((*mesh))); 5210598e1a2SLisandro Dalcin PetscFunctionReturn(0); 5220598e1a2SLisandro Dalcin } 5230598e1a2SLisandro Dalcin 5240598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh) 525de78e4feSLisandro Dalcin { 526698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 527698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 528de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 5297d5b9244SMatthew G. Knepley int n, t, num, nid, snum; 5300598e1a2SLisandro Dalcin GmshNodes *nodes; 531de78e4feSLisandro Dalcin 532de78e4feSLisandro Dalcin PetscFunctionBegin; 5339566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 534698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 53508401ef6SPierre Jolivet PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 5369566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(num, &nodes)); 5370598e1a2SLisandro Dalcin mesh->numNodes = num; 5380598e1a2SLisandro Dalcin mesh->nodelist = nodes; 5390598e1a2SLisandro Dalcin for (n = 0; n < num; ++n) { 5400598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 5419566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM)); 5429566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE)); 5439566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 5449566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 5450598e1a2SLisandro Dalcin nodes->id[n] = nid; 5467d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n*GMSH_MAX_TAGS+t] = -1; 54711c56cb0SLisandro Dalcin } 54811c56cb0SLisandro Dalcin PetscFunctionReturn(0); 54911c56cb0SLisandro Dalcin } 55011c56cb0SLisandro Dalcin 551de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 552de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 553de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 554de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 5550598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v22(GmshFile* gmsh, GmshMesh *mesh) 55611c56cb0SLisandro Dalcin { 557698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 558698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 559698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 560de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 5610598e1a2SLisandro Dalcin int i, c, p, num, ibuf[1+4+1000], snum; 5620598e1a2SLisandro Dalcin int cellType, numElem, numVerts, numNodes, numTags; 563df032b82SMatthew G. Knepley GmshElement *elements; 5640598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 565df032b82SMatthew G. Knepley 566df032b82SMatthew G. Knepley PetscFunctionBegin; 5679566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 568698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 56908401ef6SPierre Jolivet PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 5709566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(num, &elements)); 5710598e1a2SLisandro Dalcin mesh->numElems = num; 5720598e1a2SLisandro Dalcin mesh->elements = elements; 573698a79b9SLisandro Dalcin for (c = 0; c < num;) { 5749566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM)); 5759566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3)); 5760598e1a2SLisandro Dalcin 5770598e1a2SLisandro Dalcin cellType = binary ? ibuf[0] : ibuf[1]; 5780598e1a2SLisandro Dalcin numElem = binary ? ibuf[1] : 1; 579df032b82SMatthew G. Knepley numTags = ibuf[2]; 5800598e1a2SLisandro Dalcin 5819566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 5820598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 5830598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 5840598e1a2SLisandro Dalcin 585df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 5860598e1a2SLisandro Dalcin GmshElement *element = elements + c; 5870598e1a2SLisandro Dalcin const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off; 5880598e1a2SLisandro Dalcin const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1; 5899566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf+off, nint, NULL, PETSC_ENUM)); 5909566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf+off, PETSC_ENUM, nint)); 5910598e1a2SLisandro Dalcin element->id = id[0]; 5920598e1a2SLisandro Dalcin element->dim = GmshCellMap[cellType].dim; 5930598e1a2SLisandro Dalcin element->cellType = cellType; 5940598e1a2SLisandro Dalcin element->numVerts = numVerts; 5950598e1a2SLisandro Dalcin element->numNodes = numNodes; 5967d5b9244SMatthew G. Knepley element->numTags = PetscMin(numTags, GMSH_MAX_TAGS); 5979566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 5980598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 5990598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = tags[p]; 600df032b82SMatthew G. Knepley } 601df032b82SMatthew G. Knepley } 602df032b82SMatthew G. Knepley PetscFunctionReturn(0); 603df032b82SMatthew G. Knepley } 6047d282ae0SMichael Lange 605de78e4feSLisandro Dalcin /* 606de78e4feSLisandro Dalcin $Entities 607de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 608de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 609de78e4feSLisandro Dalcin // points 610de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 611de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 612de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 613de78e4feSLisandro Dalcin ... 614de78e4feSLisandro Dalcin // curves 615de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 616de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 617de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 618de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 619de78e4feSLisandro Dalcin ... 620de78e4feSLisandro Dalcin // surfaces 621de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 622de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 623de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 624de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 625de78e4feSLisandro Dalcin ... 626de78e4feSLisandro Dalcin // volumes 627de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 628de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 629de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 630de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 631de78e4feSLisandro Dalcin ... 632de78e4feSLisandro Dalcin $EndEntities 633de78e4feSLisandro Dalcin */ 6340598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh) 635de78e4feSLisandro Dalcin { 636698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 637698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 638698a79b9SLisandro Dalcin long index, num, lbuf[4]; 639730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 640698a79b9SLisandro Dalcin PetscInt count[4], i; 641a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 642de78e4feSLisandro Dalcin 643de78e4feSLisandro Dalcin PetscFunctionBegin; 6449566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG)); 6459566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(lbuf, PETSC_LONG, 4)); 646698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 6479566063dSJacob Faibussowitsch PetscCall(GmshEntitiesCreate(count, &mesh->entities)); 648de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 649730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 6509566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM)); 6519566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&eid, PETSC_ENUM, 1)); 6529566063dSJacob Faibussowitsch PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity)); 6539566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE)); 6549566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6)); 6559566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG)); 6569566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1)); 6579566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf)); 6589566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM)); 6599566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num)); 6607d5b9244SMatthew G. Knepley entity->numTags = numTags = (int) PetscMin(num, GMSH_MAX_TAGS); 661de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 662de78e4feSLisandro Dalcin if (dim == 0) continue; 6639566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG)); 6649566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&num, PETSC_LONG, 1)); 6659566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, num, sizeof(int), &ibuf)); 6669566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM)); 6679566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, num)); 668de78e4feSLisandro Dalcin } 669de78e4feSLisandro Dalcin } 670de78e4feSLisandro Dalcin PetscFunctionReturn(0); 671de78e4feSLisandro Dalcin } 672de78e4feSLisandro Dalcin 673de78e4feSLisandro Dalcin /* 674de78e4feSLisandro Dalcin $Nodes 675de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 676de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 677de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 678de78e4feSLisandro Dalcin ... 679de78e4feSLisandro Dalcin ... 680de78e4feSLisandro Dalcin $EndNodes 681de78e4feSLisandro Dalcin */ 6820598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh) 683de78e4feSLisandro Dalcin { 684698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 685698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 6867d5b9244SMatthew G. Knepley long block, node, n, t, numEntityBlocks, numTotalNodes, numNodes; 687de78e4feSLisandro Dalcin int info[3], nid; 6880598e1a2SLisandro Dalcin GmshNodes *nodes; 689de78e4feSLisandro Dalcin 690de78e4feSLisandro Dalcin PetscFunctionBegin; 6919566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG)); 6929566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1)); 6939566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG)); 6949566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numTotalNodes, PETSC_LONG, 1)); 6959566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(numTotalNodes, &nodes)); 6960598e1a2SLisandro Dalcin mesh->numNodes = numTotalNodes; 6970598e1a2SLisandro Dalcin mesh->nodelist = nodes; 6980598e1a2SLisandro Dalcin for (n = 0, block = 0; block < numEntityBlocks; ++block) { 6999566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM)); 7009566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG)); 7019566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numNodes, PETSC_LONG, 1)); 702698a79b9SLisandro Dalcin if (gmsh->binary) { 703277f51e8SBarry Smith size_t nbytes = sizeof(int) + 3*sizeof(double); 704277f51e8SBarry Smith char *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */ 7059566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numNodes, nbytes, &cbuf)); 7069566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR)); 7070598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 708de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 7090598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 7109566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cnid, PETSC_ENUM, 1)); 7119566063dSJacob Faibussowitsch if (!PetscBinaryBigEndian()) PetscCall(PetscByteSwap(cxyz, PETSC_DOUBLE, 3)); 7129566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&nid, cnid, sizeof(int))); 7139566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(xyz, cxyz, 3*sizeof(double))); 7149566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 7159566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 7160598e1a2SLisandro Dalcin nodes->id[n] = nid; 7177d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n*GMSH_MAX_TAGS+t] = -1; 718de78e4feSLisandro Dalcin } 719de78e4feSLisandro Dalcin } else { 7200598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 7210598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 7229566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM)); 7239566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE)); 7249566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nid, PETSC_ENUM, 1)); 7259566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(xyz, PETSC_DOUBLE, 3)); 7260598e1a2SLisandro Dalcin nodes->id[n] = nid; 7277d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) nodes->tag[n*GMSH_MAX_TAGS+t] = -1; 728de78e4feSLisandro Dalcin } 729de78e4feSLisandro Dalcin } 730de78e4feSLisandro Dalcin } 731de78e4feSLisandro Dalcin PetscFunctionReturn(0); 732de78e4feSLisandro Dalcin } 733de78e4feSLisandro Dalcin 734de78e4feSLisandro Dalcin /* 735de78e4feSLisandro Dalcin $Elements 736de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 737de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 738de78e4feSLisandro Dalcin tag(int) numVert[...](int) 739de78e4feSLisandro Dalcin ... 740de78e4feSLisandro Dalcin ... 741de78e4feSLisandro Dalcin $EndElements 742de78e4feSLisandro Dalcin */ 7430598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh) 744de78e4feSLisandro Dalcin { 745698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 746698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 747de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 748de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 7490598e1a2SLisandro Dalcin int eid, dim, cellType, numVerts, numNodes, numTags; 750a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 751de78e4feSLisandro Dalcin GmshElement *elements; 7520598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 753de78e4feSLisandro Dalcin 754de78e4feSLisandro Dalcin PetscFunctionBegin; 7559566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG)); 7569566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1)); 7579566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG)); 7589566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numTotalElements, PETSC_LONG, 1)); 7599566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(numTotalElements, &elements)); 7600598e1a2SLisandro Dalcin mesh->numElems = numTotalElements; 7610598e1a2SLisandro Dalcin mesh->elements = elements; 762de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 7639566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM)); 7649566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(info, PETSC_ENUM, 3)); 765de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 7669566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 7679566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 7680598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 7690598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 770730356f6SLisandro Dalcin numTags = entity->numTags; 7719566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG)); 7729566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numElements, PETSC_LONG, 1)); 7739566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf)); 7749566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM)); 7759566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements)); 776de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 777de78e4feSLisandro Dalcin GmshElement *element = elements + c; 7780598e1a2SLisandro Dalcin const int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 7790598e1a2SLisandro Dalcin element->id = id[0]; 780de78e4feSLisandro Dalcin element->dim = dim; 781de78e4feSLisandro Dalcin element->cellType = cellType; 7820598e1a2SLisandro Dalcin element->numVerts = numVerts; 783de78e4feSLisandro Dalcin element->numNodes = numNodes; 784de78e4feSLisandro Dalcin element->numTags = numTags; 7859566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 7860598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 7870598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 788de78e4feSLisandro Dalcin } 789de78e4feSLisandro Dalcin } 790698a79b9SLisandro Dalcin PetscFunctionReturn(0); 791698a79b9SLisandro Dalcin } 792698a79b9SLisandro Dalcin 7930598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[]) 794698a79b9SLisandro Dalcin { 795698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 796698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 797698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 798698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 799698a79b9SLisandro Dalcin int numPeriodic, snum, i; 800698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 8010598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 802698a79b9SLisandro Dalcin 803698a79b9SLisandro Dalcin PetscFunctionBegin; 804698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 8059566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 806698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 80708401ef6SPierre Jolivet PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 808698a79b9SLisandro Dalcin } else { 8099566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM)); 8109566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&numPeriodic, PETSC_ENUM, 1)); 811698a79b9SLisandro Dalcin } 812698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 8139dddd249SSatish Balay int ibuf[3], correspondingDim = -1, correspondingTag = -1, primaryTag = -1, correspondingNode, primaryNode; 814698a79b9SLisandro Dalcin long j, nNodes; 815698a79b9SLisandro Dalcin double affine[16]; 816698a79b9SLisandro Dalcin 817698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 8189566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING)); 8199dddd249SSatish Balay snum = sscanf(line, "%d %d %d", &correspondingDim, &correspondingTag, &primaryTag); 82008401ef6SPierre Jolivet PetscCheck(snum == 3,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 821698a79b9SLisandro Dalcin } else { 8229566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM)); 8239566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 3)); 8249dddd249SSatish Balay correspondingDim = ibuf[0]; correspondingTag = ibuf[1]; primaryTag = ibuf[2]; 825698a79b9SLisandro Dalcin } 8269dddd249SSatish Balay (void)correspondingDim; (void)correspondingTag; (void)primaryTag; /* unused */ 827698a79b9SLisandro Dalcin 828698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 8299566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 830698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 8314c500f23SPierre Jolivet if (snum != 1) { /* discard transformation and try again */ 8329566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING)); 8339566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 834698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 83508401ef6SPierre Jolivet PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 836698a79b9SLisandro Dalcin } 837698a79b9SLisandro Dalcin } else { 8389566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG)); 8399566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1)); 8404c500f23SPierre Jolivet if (nNodes == -1) { /* discard transformation and try again */ 8419566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE)); 8429566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG)); 8439566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&nNodes, PETSC_LONG, 1)); 844698a79b9SLisandro Dalcin } 845698a79b9SLisandro Dalcin } 846698a79b9SLisandro Dalcin 847698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 848698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 8499566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 8509dddd249SSatish Balay snum = sscanf(line, "%d %d", &correspondingNode, &primaryNode); 85108401ef6SPierre Jolivet PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 852698a79b9SLisandro Dalcin } else { 8539566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM)); 8549566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, 2)); 8559dddd249SSatish Balay correspondingNode = ibuf[0]; primaryNode = ibuf[1]; 856698a79b9SLisandro Dalcin } 8579dddd249SSatish Balay correspondingNode = (int) nodeMap[correspondingNode]; 8589dddd249SSatish Balay primaryNode = (int) nodeMap[primaryNode]; 8599dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 860698a79b9SLisandro Dalcin } 861698a79b9SLisandro Dalcin } 862698a79b9SLisandro Dalcin PetscFunctionReturn(0); 863698a79b9SLisandro Dalcin } 864698a79b9SLisandro Dalcin 8650598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 866698a79b9SLisandro Dalcin $Entities 867698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 868698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 869698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 870698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 871698a79b9SLisandro Dalcin ... 872698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 873698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 874698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 875698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 876698a79b9SLisandro Dalcin ... 877698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 878698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 879698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 880698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 881698a79b9SLisandro Dalcin ... 882698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 883698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 884698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 885698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 886698a79b9SLisandro Dalcin ... 887698a79b9SLisandro Dalcin $EndEntities 888698a79b9SLisandro Dalcin */ 8890598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh) 890698a79b9SLisandro Dalcin { 891698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 892698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 893698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 894698a79b9SLisandro Dalcin 895698a79b9SLisandro Dalcin PetscFunctionBegin; 8969566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, count, 4)); 8979566063dSJacob Faibussowitsch PetscCall(GmshEntitiesCreate(count, &mesh->entities)); 898698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 899698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 9009566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, &eid, 1)); 9019566063dSJacob Faibussowitsch PetscCall(GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity)); 9029566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6)); 9039566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numTags, 1)); 9049566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags)); 9059566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, tags, numTags)); 9069c0edc59SMatthew 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); 9077d5b9244SMatthew G. Knepley entity->numTags = numTags; 908698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 909698a79b9SLisandro Dalcin if (dim == 0) continue; 9109566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numTags, 1)); 9119566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numTags, sizeof(int), &tags)); 9129566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, tags, numTags)); 91381a1af93SMatthew G. Knepley /* Currently, we do not save the ids for the bounding entities */ 914698a79b9SLisandro Dalcin } 915698a79b9SLisandro Dalcin } 916698a79b9SLisandro Dalcin PetscFunctionReturn(0); 917698a79b9SLisandro Dalcin } 918698a79b9SLisandro Dalcin 91933a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 920698a79b9SLisandro Dalcin $Nodes 921698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 922698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 923698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 924698a79b9SLisandro Dalcin nodeTag(size_t) 925698a79b9SLisandro Dalcin ... 926698a79b9SLisandro Dalcin x(double) y(double) z(double) 927698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 928698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 929698a79b9SLisandro Dalcin ... 930698a79b9SLisandro Dalcin ... 931698a79b9SLisandro Dalcin $EndNodes 932698a79b9SLisandro Dalcin */ 9330598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh) 934698a79b9SLisandro Dalcin { 93581a1af93SMatthew G. Knepley int info[3], dim, eid, parametric; 9367d5b9244SMatthew G. Knepley PetscInt sizes[4], numEntityBlocks, numTags, t, numNodes, numNodesBlock = 0, block, node, n; 93781a1af93SMatthew G. Knepley GmshEntity *entity = NULL; 9380598e1a2SLisandro Dalcin GmshNodes *nodes; 939698a79b9SLisandro Dalcin 940698a79b9SLisandro Dalcin PetscFunctionBegin; 9419566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, sizes, 4)); 9420598e1a2SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 9439566063dSJacob Faibussowitsch PetscCall(GmshNodesCreate(numNodes, &nodes)); 9440598e1a2SLisandro Dalcin mesh->numNodes = numNodes; 9450598e1a2SLisandro Dalcin mesh->nodelist = nodes; 946698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 9479566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 94881a1af93SMatthew G. Knepley dim = info[0]; eid = info[1]; parametric = info[2]; 9499566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 9507d5b9244SMatthew G. Knepley numTags = entity->numTags; 95181a1af93SMatthew G. Knepley PetscCheck(!parametric, PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 9529566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numNodesBlock, 1)); 9539566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, nodes->id+node, numNodesBlock)); 9549566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, nodes->xyz+node*3, numNodesBlock*3)); 9557d5b9244SMatthew G. Knepley for (n = 0; n < numNodesBlock; ++n) { 9567d5b9244SMatthew G. Knepley PetscInt *tags = &nodes->tag[node*GMSH_MAX_TAGS]; 9577d5b9244SMatthew G. Knepley 9587d5b9244SMatthew G. Knepley for (t = 0; t < numTags; ++t) tags[n*GMSH_MAX_TAGS+t] = entity->tags[t]; 9597d5b9244SMatthew G. Knepley for (t = numTags; t < GMSH_MAX_TAGS; ++t) tags[n*GMSH_MAX_TAGS+t] = -1; 9607d5b9244SMatthew G. Knepley } 961698a79b9SLisandro Dalcin } 9620598e1a2SLisandro Dalcin gmsh->nodeStart = sizes[2]; 9630598e1a2SLisandro Dalcin gmsh->nodeEnd = sizes[3]+1; 964698a79b9SLisandro Dalcin PetscFunctionReturn(0); 965698a79b9SLisandro Dalcin } 966698a79b9SLisandro Dalcin 96733a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 968698a79b9SLisandro Dalcin $Elements 969698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 970698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 971698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 972698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 973698a79b9SLisandro Dalcin ... 974698a79b9SLisandro Dalcin ... 975698a79b9SLisandro Dalcin $EndElements 976698a79b9SLisandro Dalcin */ 9770598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh) 978698a79b9SLisandro Dalcin { 9790598e1a2SLisandro Dalcin int info[3], eid, dim, cellType; 9800598e1a2SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p; 981698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 982698a79b9SLisandro Dalcin GmshElement *elements; 9830598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 984698a79b9SLisandro Dalcin 985698a79b9SLisandro Dalcin PetscFunctionBegin; 9869566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, sizes, 4)); 987698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 9889566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(numElements, &elements)); 9890598e1a2SLisandro Dalcin mesh->numElems = numElements; 9900598e1a2SLisandro Dalcin mesh->elements = elements; 991698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 9929566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 993698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 9949566063dSJacob Faibussowitsch PetscCall(GmshEntitiesGet(mesh->entities, dim, eid, &entity)); 9959566063dSJacob Faibussowitsch PetscCall(GmshCellTypeCheck(cellType)); 9960598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 9970598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 9987d5b9244SMatthew G. Knepley numTags = entity->numTags; 9999566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numBlockElements, 1)); 10009566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf)); 10019566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements)); 1002698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 1003698a79b9SLisandro Dalcin GmshElement *element = elements + c; 10040598e1a2SLisandro Dalcin const PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 1005698a79b9SLisandro Dalcin element->id = id[0]; 1006698a79b9SLisandro Dalcin element->dim = dim; 1007698a79b9SLisandro Dalcin element->cellType = cellType; 10080598e1a2SLisandro Dalcin element->numVerts = numVerts; 1009698a79b9SLisandro Dalcin element->numNodes = numNodes; 1010698a79b9SLisandro Dalcin element->numTags = numTags; 10119566063dSJacob Faibussowitsch PetscCall(PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes)); 10120598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 10130598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 1014698a79b9SLisandro Dalcin } 1015698a79b9SLisandro Dalcin } 1016698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1017698a79b9SLisandro Dalcin } 1018698a79b9SLisandro Dalcin 10190598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 1020698a79b9SLisandro Dalcin $Periodic 1021698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 10229dddd249SSatish Balay entityDim(int) entityTag(int) entityTagPrimary(int) 1023698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 1024698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 10259dddd249SSatish Balay nodeTag(size_t) nodeTagPrimary(size_t) 1026698a79b9SLisandro Dalcin ... 1027698a79b9SLisandro Dalcin ... 1028698a79b9SLisandro Dalcin $EndPeriodic 1029698a79b9SLisandro Dalcin */ 10300598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[]) 1031698a79b9SLisandro Dalcin { 1032698a79b9SLisandro Dalcin int info[3]; 1033698a79b9SLisandro Dalcin double dbuf[16]; 10340598e1a2SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 10350598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 1036698a79b9SLisandro Dalcin 1037698a79b9SLisandro Dalcin PetscFunctionBegin; 10389566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numPeriodicLinks, 1)); 1039698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 10409566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, info, 3)); 10419566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numAffine, 1)); 10429566063dSJacob Faibussowitsch PetscCall(GmshReadDouble(gmsh, dbuf, numAffine)); 10439566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, &numCorrespondingNodes, 1)); 10449566063dSJacob Faibussowitsch PetscCall(GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags)); 10459566063dSJacob Faibussowitsch PetscCall(GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2)); 1046698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 10479dddd249SSatish Balay PetscInt correspondingNode = nodeMap[nodeTags[node*2+0]]; 10489dddd249SSatish Balay PetscInt primaryNode = nodeMap[nodeTags[node*2+1]]; 10499dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 1050698a79b9SLisandro Dalcin } 1051698a79b9SLisandro Dalcin } 1052698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1053698a79b9SLisandro Dalcin } 1054698a79b9SLisandro Dalcin 10550598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 1056d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 1057d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 1058d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 1059d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 1060d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 1061d6689ca9SLisandro Dalcin $EndMeshFormat 1062d6689ca9SLisandro Dalcin */ 10630598e1a2SLisandro Dalcin static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh) 1064698a79b9SLisandro Dalcin { 1065698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1066698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 1067698a79b9SLisandro Dalcin float version; 1068698a79b9SLisandro Dalcin 1069698a79b9SLisandro Dalcin PetscFunctionBegin; 10709566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 3)); 1071698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 1072a8d4e440SJunchao Zhang fileFormat = (int)roundf(version*10); 107308401ef6SPierre Jolivet PetscCheck(snum == 3,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 1074a8d4e440SJunchao Zhang PetscCheck(fileFormat >= 22,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 10751dca8a05SBarry Smith PetscCheck((int)version != 3,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 1076a8d4e440SJunchao Zhang PetscCheck(fileFormat <= 41,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 107708401ef6SPierre Jolivet PetscCheck(!gmsh->binary || fileType,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII"); 107808401ef6SPierre Jolivet PetscCheck(gmsh->binary || !fileType,PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary"); 10791dca8a05SBarry 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); 10801dca8a05SBarry 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); 1081698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 1082698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 1083698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 1084698a79b9SLisandro Dalcin if (gmsh->binary) { 10859566063dSJacob Faibussowitsch PetscCall(GmshReadInt(gmsh, &checkEndian, 1)); 1086698a79b9SLisandro Dalcin if (checkEndian != 1) { 10879566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(&checkEndian, PETSC_ENUM, 1)); 108808401ef6SPierre Jolivet PetscCheck(checkEndian == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line); 1089698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 1090698a79b9SLisandro Dalcin } 1091698a79b9SLisandro Dalcin } 1092698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1093698a79b9SLisandro Dalcin } 1094698a79b9SLisandro Dalcin 10951f643af3SLisandro Dalcin /* 10961f643af3SLisandro Dalcin PhysicalNames 10971f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 10981f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 10991f643af3SLisandro Dalcin ... 11001f643af3SLisandro Dalcin $EndPhysicalNames 11011f643af3SLisandro Dalcin */ 1102a45dabc8SMatthew G. Knepley static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh, GmshMesh *mesh) 1103698a79b9SLisandro Dalcin { 11045f5cd6d5SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q, *r; 1105a45dabc8SMatthew G. Knepley int snum, region, dim, tag; 1106698a79b9SLisandro Dalcin 1107698a79b9SLisandro Dalcin PetscFunctionBegin; 11089566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 1)); 1109a45dabc8SMatthew G. Knepley snum = sscanf(line, "%d", ®ion); 1110a45dabc8SMatthew G. Knepley mesh->numRegions = region; 111108401ef6SPierre Jolivet PetscCheck(snum == 1,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 1112*90d778b8SLisandro Dalcin PetscCall(PetscMalloc3(mesh->numRegions, &mesh->regionDims, mesh->numRegions, &mesh->regionTags, mesh->numRegions, &mesh->regionNames)); 1113a45dabc8SMatthew G. Knepley for (region = 0; region < mesh->numRegions; ++region) { 11149566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 2)); 11151f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 111608401ef6SPierre Jolivet PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11179566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, -(PetscInt)sizeof(line))); 11189566063dSJacob Faibussowitsch PetscCall(PetscStrchr(line, '"', &p)); 111928b400f6SJacob Faibussowitsch PetscCheck(p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11209566063dSJacob Faibussowitsch PetscCall(PetscStrrchr(line, '"', &q)); 112108401ef6SPierre Jolivet PetscCheck(q != p, PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11225f5cd6d5SMatthew G. Knepley PetscCall(PetscStrrchr(line, ':', &r)); 11235f5cd6d5SMatthew G. Knepley if (p != r) q = r; 11249566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(name, p+1, (size_t)(q-p-1))); 1125*90d778b8SLisandro Dalcin mesh->regionDims[region] = dim; 1126a45dabc8SMatthew G. Knepley mesh->regionTags[region] = tag; 11279566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &mesh->regionNames[region])); 1128698a79b9SLisandro Dalcin } 1129de78e4feSLisandro Dalcin PetscFunctionReturn(0); 1130de78e4feSLisandro Dalcin } 1131de78e4feSLisandro Dalcin 11320598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh) 113396ca5757SLisandro Dalcin { 11340598e1a2SLisandro Dalcin PetscFunctionBegin; 11350598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11369566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadEntities_v41(gmsh, mesh)); break; 11379566063dSJacob Faibussowitsch default: PetscCall(GmshReadEntities_v40(gmsh, mesh)); break; 113896ca5757SLisandro Dalcin } 11390598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11400598e1a2SLisandro Dalcin } 11410598e1a2SLisandro Dalcin 11420598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh) 11430598e1a2SLisandro Dalcin { 11440598e1a2SLisandro Dalcin PetscFunctionBegin; 11450598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11469566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadNodes_v41(gmsh, mesh)); break; 11479566063dSJacob Faibussowitsch case 40: PetscCall(GmshReadNodes_v40(gmsh, mesh)); break; 11489566063dSJacob Faibussowitsch default: PetscCall(GmshReadNodes_v22(gmsh, mesh)); break; 11490598e1a2SLisandro Dalcin } 11500598e1a2SLisandro Dalcin 11510598e1a2SLisandro Dalcin { /* Gmsh v2.2/v4.0 does not provide min/max node tags */ 11520598e1a2SLisandro Dalcin if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) { 11530598e1a2SLisandro Dalcin PetscInt tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n; 11540598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11550598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11560598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 11570598e1a2SLisandro Dalcin tagMin = PetscMin(tag, tagMin); 11580598e1a2SLisandro Dalcin tagMax = PetscMax(tag, tagMax); 11590598e1a2SLisandro Dalcin } 11600598e1a2SLisandro Dalcin gmsh->nodeStart = tagMin; 11610598e1a2SLisandro Dalcin gmsh->nodeEnd = tagMax+1; 11620598e1a2SLisandro Dalcin } 11630598e1a2SLisandro Dalcin } 11640598e1a2SLisandro Dalcin 11650598e1a2SLisandro Dalcin { /* Support for sparse node tags */ 11660598e1a2SLisandro Dalcin PetscInt n, t; 11670598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf)); 11690598e1a2SLisandro Dalcin for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT; 11700598e1a2SLisandro Dalcin gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart; 11710598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11720598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 117363a3b9bcSJacob Faibussowitsch PetscCheck(gmsh->nodeMap[tag] < 0,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %" PetscInt_FMT, tag); 11740598e1a2SLisandro Dalcin gmsh->nodeMap[tag] = n; 11750598e1a2SLisandro Dalcin } 11760598e1a2SLisandro Dalcin } 11770598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11780598e1a2SLisandro Dalcin } 11790598e1a2SLisandro Dalcin 11800598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh) 11810598e1a2SLisandro Dalcin { 11820598e1a2SLisandro Dalcin PetscFunctionBegin; 11830598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11849566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadElements_v41(gmsh, mesh)); break; 11859566063dSJacob Faibussowitsch case 40: PetscCall(GmshReadElements_v40(gmsh, mesh)); break; 11869566063dSJacob Faibussowitsch default: PetscCall(GmshReadElements_v22(gmsh, mesh)); break; 11870598e1a2SLisandro Dalcin } 11880598e1a2SLisandro Dalcin 11890598e1a2SLisandro Dalcin { /* Reorder elements by codimension and polytope type */ 11900598e1a2SLisandro Dalcin PetscInt ne = mesh->numElems; 11910598e1a2SLisandro Dalcin GmshElement *elements = mesh->elements; 1192066ea43fSLisandro Dalcin PetscInt keymap[GMSH_NUM_POLYTOPES], nk = 0; 1193066ea43fSLisandro Dalcin PetscInt offset[GMSH_NUM_POLYTOPES+1], e, k; 11940598e1a2SLisandro Dalcin 1195066ea43fSLisandro Dalcin for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT; 11969566063dSJacob Faibussowitsch PetscCall(PetscMemzero(offset,sizeof(offset))); 11970598e1a2SLisandro Dalcin 1198066ea43fSLisandro Dalcin keymap[GMSH_TET] = nk++; 1199066ea43fSLisandro Dalcin keymap[GMSH_HEX] = nk++; 1200066ea43fSLisandro Dalcin keymap[GMSH_PRI] = nk++; 1201066ea43fSLisandro Dalcin keymap[GMSH_PYR] = nk++; 1202066ea43fSLisandro Dalcin keymap[GMSH_TRI] = nk++; 1203066ea43fSLisandro Dalcin keymap[GMSH_QUA] = nk++; 1204066ea43fSLisandro Dalcin keymap[GMSH_SEG] = nk++; 1205066ea43fSLisandro Dalcin keymap[GMSH_VTX] = nk++; 12060598e1a2SLisandro Dalcin 12079566063dSJacob Faibussowitsch PetscCall(GmshElementsCreate(mesh->numElems, &mesh->elements)); 12080598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope] 12090598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) offset[1+key(e)]++; 12100598e1a2SLisandro Dalcin for (k = 1; k < nk; ++k) offset[k] += offset[k-1]; 12110598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e]; 12120598e1a2SLisandro Dalcin #undef key 12139566063dSJacob Faibussowitsch PetscCall(GmshElementsDestroy(&elements)); 1214066ea43fSLisandro Dalcin } 12150598e1a2SLisandro Dalcin 1216066ea43fSLisandro Dalcin { /* Mesh dimension and order */ 1217066ea43fSLisandro Dalcin GmshElement *elem = mesh->numElems ? mesh->elements : NULL; 1218066ea43fSLisandro Dalcin mesh->dim = elem ? GmshCellMap[elem->cellType].dim : 0; 1219066ea43fSLisandro Dalcin mesh->order = elem ? GmshCellMap[elem->cellType].order : 0; 12200598e1a2SLisandro Dalcin } 12210598e1a2SLisandro Dalcin 12220598e1a2SLisandro Dalcin { 12230598e1a2SLisandro Dalcin PetscBT vtx; 12240598e1a2SLisandro Dalcin PetscInt dim = mesh->dim, e, n, v; 12250598e1a2SLisandro Dalcin 12269566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(mesh->numNodes, &vtx)); 12270598e1a2SLisandro Dalcin 12280598e1a2SLisandro Dalcin /* Compute number of cells and set of vertices */ 12290598e1a2SLisandro Dalcin mesh->numCells = 0; 12300598e1a2SLisandro Dalcin for (e = 0; e < mesh->numElems; ++e) { 12310598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 12320598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) mesh->numCells++; 12330598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; v++) { 12349566063dSJacob Faibussowitsch PetscCall(PetscBTSet(vtx, elem->nodes[v])); 12350598e1a2SLisandro Dalcin } 12360598e1a2SLisandro Dalcin } 12370598e1a2SLisandro Dalcin 12380598e1a2SLisandro Dalcin /* Compute numbering for vertices */ 12390598e1a2SLisandro Dalcin mesh->numVerts = 0; 12409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mesh->numNodes, &mesh->vertexMap)); 12410598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12420598e1a2SLisandro Dalcin mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT; 12430598e1a2SLisandro Dalcin 12449566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&vtx)); 12450598e1a2SLisandro Dalcin } 12460598e1a2SLisandro Dalcin PetscFunctionReturn(0); 12470598e1a2SLisandro Dalcin } 12480598e1a2SLisandro Dalcin 12490598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh) 12500598e1a2SLisandro Dalcin { 12510598e1a2SLisandro Dalcin PetscInt n; 12520598e1a2SLisandro Dalcin 12530598e1a2SLisandro Dalcin PetscFunctionBegin; 12549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(mesh->numNodes, &mesh->periodMap)); 12550598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n; 12560598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 12579566063dSJacob Faibussowitsch case 41: PetscCall(GmshReadPeriodic_v41(gmsh, mesh->periodMap)); break; 12589566063dSJacob Faibussowitsch default: PetscCall(GmshReadPeriodic_v40(gmsh, mesh->periodMap)); break; 12590598e1a2SLisandro Dalcin } 12600598e1a2SLisandro Dalcin 12619dddd249SSatish Balay /* Find canonical primary nodes */ 12620598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12630598e1a2SLisandro Dalcin while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]]) 12640598e1a2SLisandro Dalcin mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]]; 12650598e1a2SLisandro Dalcin 12669dddd249SSatish Balay /* Renumber vertices (filter out correspondings) */ 12670598e1a2SLisandro Dalcin mesh->numVerts = 0; 12680598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12690598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12709dddd249SSatish Balay if (mesh->periodMap[n] == n) /* is primary */ 12710598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->numVerts++; 12720598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12730598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12749dddd249SSatish Balay if (mesh->periodMap[n] != n) /* is corresponding */ 12750598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]]; 12760598e1a2SLisandro Dalcin PetscFunctionReturn(0); 12770598e1a2SLisandro Dalcin } 12780598e1a2SLisandro Dalcin 1279066ea43fSLisandro Dalcin #define DM_POLYTOPE_VERTEX DM_POLYTOPE_POINT 1280066ea43fSLisandro Dalcin static const DMPolytopeType DMPolytopeMap[] = { 1281066ea43fSLisandro Dalcin /* GMSH_VTX */ DM_POLYTOPE_VERTEX, 1282066ea43fSLisandro Dalcin /* GMSH_SEG */ DM_POLYTOPE_SEGMENT, 1283066ea43fSLisandro Dalcin /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE, 1284066ea43fSLisandro Dalcin /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL, 1285066ea43fSLisandro Dalcin /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON, 1286066ea43fSLisandro Dalcin /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON, 1287066ea43fSLisandro Dalcin /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM, 1288066ea43fSLisandro Dalcin /* GMSH_PYR */ DM_POLYTOPE_PYRAMID, 1289066ea43fSLisandro Dalcin DM_POLYTOPE_UNKNOWN 1290066ea43fSLisandro Dalcin }; 12910598e1a2SLisandro Dalcin 12929fbee547SJacob Faibussowitsch static inline DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType) 12930598e1a2SLisandro Dalcin { 1294066ea43fSLisandro Dalcin return DMPolytopeMap[GmshCellMap[cellType].polytope]; 1295066ea43fSLisandro Dalcin } 1296066ea43fSLisandro Dalcin 1297066ea43fSLisandro Dalcin static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem) 1298066ea43fSLisandro Dalcin { 1299066ea43fSLisandro Dalcin DM K; 1300066ea43fSLisandro Dalcin PetscSpace P; 1301066ea43fSLisandro Dalcin PetscDualSpace Q; 1302066ea43fSLisandro Dalcin PetscQuadrature q, fq; 1303066ea43fSLisandro Dalcin PetscBool isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 1304066ea43fSLisandro Dalcin PetscBool endpoint = PETSC_TRUE; 1305066ea43fSLisandro Dalcin char name[32]; 1306066ea43fSLisandro Dalcin 1307066ea43fSLisandro Dalcin PetscFunctionBegin; 1308066ea43fSLisandro Dalcin /* Create space */ 13099566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(comm, &P)); 13109566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL)); 13119566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialSetTensor(P, isTensor)); 13129566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(P, Nc)); 13139566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(P, dim)); 13149566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(P, k, PETSC_DETERMINE)); 1315066ea43fSLisandro Dalcin if (prefix) { 13169566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, prefix)); 13179566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetFromOptions(P)); 13189566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) P, NULL)); 13199566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(P, &k, NULL)); 1320066ea43fSLisandro Dalcin } 13219566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetUp(P)); 1322066ea43fSLisandro Dalcin /* Create dual space */ 13239566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceCreate(comm, &Q)); 13249566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE)); 13259566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetTensor(Q, isTensor)); 13269566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetContinuity(Q, continuity)); 13279566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0)); 13289566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetNumComponents(Q, Nc)); 13299566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetOrder(Q, k)); 13309566063dSJacob Faibussowitsch PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, isSimplex), &K)); 13319566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetDM(Q, K)); 13329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&K)); 1333066ea43fSLisandro Dalcin if (prefix) { 13349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, prefix)); 13359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetFromOptions(Q)); 13369566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) Q, NULL)); 1337066ea43fSLisandro Dalcin } 13389566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceSetUp(Q)); 1339066ea43fSLisandro Dalcin /* Create quadrature */ 1340066ea43fSLisandro Dalcin if (isSimplex) { 13419566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim, 1, k+1, -1, +1, &q)); 13429566063dSJacob Faibussowitsch PetscCall(PetscDTStroudConicalQuadrature(dim-1, 1, k+1, -1, +1, &fq)); 1343066ea43fSLisandro Dalcin } else { 13449566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, k+1, -1, +1, &q)); 13459566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim-1, 1, k+1, -1, +1, &fq)); 1346066ea43fSLisandro Dalcin } 1347066ea43fSLisandro Dalcin /* Create finite element */ 13489566063dSJacob Faibussowitsch PetscCall(PetscFECreate(comm, fem)); 13499566063dSJacob Faibussowitsch PetscCall(PetscFESetType(*fem, PETSCFEBASIC)); 13509566063dSJacob Faibussowitsch PetscCall(PetscFESetNumComponents(*fem, Nc)); 13519566063dSJacob Faibussowitsch PetscCall(PetscFESetBasisSpace(*fem, P)); 13529566063dSJacob Faibussowitsch PetscCall(PetscFESetDualSpace(*fem, Q)); 13539566063dSJacob Faibussowitsch PetscCall(PetscFESetQuadrature(*fem, q)); 13549566063dSJacob Faibussowitsch PetscCall(PetscFESetFaceQuadrature(*fem, fq)); 1355066ea43fSLisandro Dalcin if (prefix) { 13569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix)); 13579566063dSJacob Faibussowitsch PetscCall(PetscFESetFromOptions(*fem)); 13589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject) *fem, NULL)); 1359066ea43fSLisandro Dalcin } 13609566063dSJacob Faibussowitsch PetscCall(PetscFESetUp(*fem)); 1361066ea43fSLisandro Dalcin /* Cleanup */ 13629566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&P)); 13639566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&Q)); 13649566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&q)); 13659566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&fq)); 1366066ea43fSLisandro Dalcin /* Set finite element name */ 136763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(name, sizeof(name), "%s%" PetscInt_FMT, isSimplex? "P" : "Q", k)); 13689566063dSJacob Faibussowitsch PetscCall(PetscFESetName(*fem, name)); 1369066ea43fSLisandro Dalcin PetscFunctionReturn(0); 137096ca5757SLisandro Dalcin } 137196ca5757SLisandro Dalcin 1372d6689ca9SLisandro Dalcin /*@C 1373d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 1374d6689ca9SLisandro Dalcin 1375d6689ca9SLisandro Dalcin + comm - The MPI communicator 1376d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 1377d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 1378d6689ca9SLisandro Dalcin 1379d6689ca9SLisandro Dalcin Output Parameter: 1380d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 1381d6689ca9SLisandro Dalcin 1382d6689ca9SLisandro Dalcin Level: beginner 1383d6689ca9SLisandro Dalcin 1384db781477SPatrick Sanan .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 1385d6689ca9SLisandro Dalcin @*/ 1386d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1387d6689ca9SLisandro Dalcin { 1388d6689ca9SLisandro Dalcin PetscViewer viewer; 1389d6689ca9SLisandro Dalcin PetscMPIInt rank; 1390d6689ca9SLisandro Dalcin int fileType; 1391d6689ca9SLisandro Dalcin PetscViewerType vtype; 1392d6689ca9SLisandro Dalcin 1393d6689ca9SLisandro Dalcin PetscFunctionBegin; 13949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1395d6689ca9SLisandro Dalcin 1396d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 1397dd400576SPatrick Sanan if (rank == 0) { 13980598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1399d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1400d6689ca9SLisandro Dalcin int snum; 1401d6689ca9SLisandro Dalcin float version; 1402a8d4e440SJunchao Zhang int fileFormat; 1403d6689ca9SLisandro Dalcin 14049566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gmsh,1)); 14059566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer)); 14069566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII)); 14079566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ)); 14089566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(gmsh->viewer, filename)); 1409d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 14109566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 14119566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$MeshFormat", line)); 14129566063dSJacob Faibussowitsch PetscCall(GmshReadString(gmsh, line, 2)); 1413d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 1414a8d4e440SJunchao Zhang fileFormat = (int)roundf(version*10); 141508401ef6SPierre Jolivet PetscCheck(snum == 2,PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 1416a8d4e440SJunchao Zhang PetscCheck(fileFormat >= 22,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 14171dca8a05SBarry Smith PetscCheck((int)version != 3,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 1418a8d4e440SJunchao Zhang PetscCheck(fileFormat <= 41,PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 14199566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&gmsh->viewer)); 1420d6689ca9SLisandro Dalcin } 14219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&fileType, 1, MPI_INT, 0, comm)); 1422d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 1423d6689ca9SLisandro Dalcin 1424d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 14259566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 14269566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, vtype)); 14279566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 14289566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 14299566063dSJacob Faibussowitsch PetscCall(DMPlexCreateGmsh(comm, viewer, interpolate, dm)); 14309566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1431d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 1432d6689ca9SLisandro Dalcin } 1433d6689ca9SLisandro Dalcin 1434331830f3SMatthew G. Knepley /*@ 14357d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 1436331830f3SMatthew G. Knepley 1437d083f849SBarry Smith Collective 1438331830f3SMatthew G. Knepley 1439331830f3SMatthew G. Knepley Input Parameters: 1440331830f3SMatthew G. Knepley + comm - The MPI communicator 1441331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1442331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1443331830f3SMatthew G. Knepley 1444331830f3SMatthew G. Knepley Output Parameter: 1445331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1446331830f3SMatthew G. Knepley 1447698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1448331830f3SMatthew G. Knepley 1449331830f3SMatthew G. Knepley Level: beginner 1450331830f3SMatthew G. Knepley 1451db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()` 1452331830f3SMatthew G. Knepley @*/ 1453331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1454331830f3SMatthew G. Knepley { 14550598e1a2SLisandro Dalcin GmshMesh *mesh = NULL; 145611c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 14570598e1a2SLisandro Dalcin PetscBT periodicVerts = NULL; 14580598e1a2SLisandro Dalcin PetscBT periodicCells = NULL; 14596858538eSMatthew G. Knepley DM cdm, cdmCell = NULL; 14606858538eSMatthew G. Knepley PetscSection cs, csCell = NULL; 14616858538eSMatthew G. Knepley Vec coordinates, coordinatesCell; 1462a45dabc8SMatthew G. Knepley DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL, marker = NULL, *regionSets; 1463066ea43fSLisandro Dalcin PetscInt dim = 0, coordDim = -1, order = 0; 14640598e1a2SLisandro Dalcin PetscInt numNodes = 0, numElems = 0, numVerts = 0, numCells = 0; 1465066ea43fSLisandro Dalcin PetscInt cell, cone[8], e, n, v, d; 146681a1af93SMatthew G. Knepley PetscBool binary, usemarker = PETSC_FALSE, useregions = PETSC_FALSE, markvertices = PETSC_FALSE; 14670598e1a2SLisandro Dalcin PetscBool hybrid = interpolate, periodic = PETSC_TRUE; 1468066ea43fSLisandro Dalcin PetscBool highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE; 1469066ea43fSLisandro Dalcin PetscBool isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE; 14700598e1a2SLisandro Dalcin PetscMPIInt rank; 1471331830f3SMatthew G. Knepley 1472331830f3SMatthew G. Knepley PetscFunctionBegin; 14739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1474d0609cedSBarry Smith PetscObjectOptionsBegin((PetscObject)viewer); 1475d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"DMPlex Gmsh options"); 14769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", hybrid, &hybrid, NULL)); 14779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL)); 14789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_highorder","Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet)); 14799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL)); 14809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL)); 14819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_use_regions", "Generate labels with region names", "DMPlexCreateGmsh", useregions, &useregions, NULL)); 14829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-dm_plex_gmsh_mark_vertices", "Add vertices to generated labels", "DMPlexCreateGmsh", markvertices, &markvertices, NULL)); 14839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE)); 1484d0609cedSBarry Smith PetscOptionsHeadEnd(); 1485d0609cedSBarry Smith PetscOptionsEnd(); 14860598e1a2SLisandro Dalcin 14879566063dSJacob Faibussowitsch PetscCall(GmshCellInfoSetUp()); 148811c56cb0SLisandro Dalcin 14899566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 14909566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 14919566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL)); 149211c56cb0SLisandro Dalcin 14939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary)); 149411c56cb0SLisandro Dalcin 149511c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 14963b46f529SLisandro Dalcin if (binary) { 149711c56cb0SLisandro Dalcin parentviewer = viewer; 14989566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer)); 149911c56cb0SLisandro Dalcin } 150011c56cb0SLisandro Dalcin 1501dd400576SPatrick Sanan if (rank == 0) { 15020598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1503698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1504698a79b9SLisandro Dalcin PetscBool match; 1505331830f3SMatthew G. Knepley 15069566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(gmsh,1)); 1507698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1508698a79b9SLisandro Dalcin gmsh->binary = binary; 1509698a79b9SLisandro Dalcin 15109566063dSJacob Faibussowitsch PetscCall(GmshMeshCreate(&mesh)); 15110598e1a2SLisandro Dalcin 1512698a79b9SLisandro Dalcin /* Read mesh format */ 15139566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 15149566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$MeshFormat", line)); 15159566063dSJacob Faibussowitsch PetscCall(GmshReadMeshFormat(gmsh)); 15169566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndMeshFormat", line)); 151711c56cb0SLisandro Dalcin 1518dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 15199566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 15209566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$PhysicalNames", line, &match)); 1521dc0ede3bSMatthew G. Knepley if (match) { 15229566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$PhysicalNames", line)); 15239566063dSJacob Faibussowitsch PetscCall(GmshReadPhysicalNames(gmsh, mesh)); 15249566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndPhysicalNames", line)); 1525698a79b9SLisandro Dalcin /* Initial read for entity section */ 15269566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 1527dc0ede3bSMatthew G. Knepley } 152811c56cb0SLisandro Dalcin 1529de78e4feSLisandro Dalcin /* Read entities */ 1530698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 15319566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Entities", line)); 15329566063dSJacob Faibussowitsch PetscCall(GmshReadEntities(gmsh, mesh)); 15339566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndEntities", line)); 1534698a79b9SLisandro Dalcin /* Initial read for nodes section */ 15359566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 1536de78e4feSLisandro Dalcin } 1537de78e4feSLisandro Dalcin 1538698a79b9SLisandro Dalcin /* Read nodes */ 15399566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Nodes", line)); 15409566063dSJacob Faibussowitsch PetscCall(GmshReadNodes(gmsh, mesh)); 15419566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndNodes", line)); 154211c56cb0SLisandro Dalcin 1543698a79b9SLisandro Dalcin /* Read elements */ 15449566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 15459566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Elements", line)); 15469566063dSJacob Faibussowitsch PetscCall(GmshReadElements(gmsh, mesh)); 15479566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndElements", line)); 1548de78e4feSLisandro Dalcin 15490598e1a2SLisandro Dalcin /* Read periodic section (OPTIONAL) */ 1550698a79b9SLisandro Dalcin if (periodic) { 15519566063dSJacob Faibussowitsch PetscCall(GmshReadSection(gmsh, line)); 15529566063dSJacob Faibussowitsch PetscCall(GmshMatch(gmsh, "$Periodic", line, &periodic)); 1553ef12879bSLisandro Dalcin } 1554ef12879bSLisandro Dalcin if (periodic) { 15559566063dSJacob Faibussowitsch PetscCall(GmshExpect(gmsh, "$Periodic", line)); 15569566063dSJacob Faibussowitsch PetscCall(GmshReadPeriodic(gmsh, mesh)); 15579566063dSJacob Faibussowitsch PetscCall(GmshReadEndSection(gmsh, "$EndPeriodic", line)); 1558698a79b9SLisandro Dalcin } 1559698a79b9SLisandro Dalcin 15609566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->wbuf)); 15619566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->sbuf)); 15629566063dSJacob Faibussowitsch PetscCall(PetscFree(gmsh->nbuf)); 15630598e1a2SLisandro Dalcin 15640598e1a2SLisandro Dalcin dim = mesh->dim; 1565066ea43fSLisandro Dalcin order = mesh->order; 15660598e1a2SLisandro Dalcin numNodes = mesh->numNodes; 15670598e1a2SLisandro Dalcin numElems = mesh->numElems; 15680598e1a2SLisandro Dalcin numVerts = mesh->numVerts; 15690598e1a2SLisandro Dalcin numCells = mesh->numCells; 1570066ea43fSLisandro Dalcin 1571066ea43fSLisandro Dalcin { 1572066ea43fSLisandro Dalcin GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL; 1573066ea43fSLisandro Dalcin GmshElement *elemB = elemA ? elemA + mesh->numCells - 1 : NULL; 1574066ea43fSLisandro Dalcin int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1; 1575066ea43fSLisandro Dalcin int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1; 1576066ea43fSLisandro Dalcin isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE; 1577066ea43fSLisandro Dalcin isHybrid = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE; 1578066ea43fSLisandro Dalcin hasTetra = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE; 1579066ea43fSLisandro Dalcin } 1580698a79b9SLisandro Dalcin } 1581698a79b9SLisandro Dalcin 1582698a79b9SLisandro Dalcin if (parentviewer) { 15839566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer)); 1584698a79b9SLisandro Dalcin } 1585698a79b9SLisandro Dalcin 1586066ea43fSLisandro Dalcin { 1587066ea43fSLisandro Dalcin int buf[6]; 1588698a79b9SLisandro Dalcin 1589066ea43fSLisandro Dalcin buf[0] = (int)dim; 1590066ea43fSLisandro Dalcin buf[1] = (int)order; 1591066ea43fSLisandro Dalcin buf[2] = periodic; 1592066ea43fSLisandro Dalcin buf[3] = isSimplex; 1593066ea43fSLisandro Dalcin buf[4] = isHybrid; 1594066ea43fSLisandro Dalcin buf[5] = hasTetra; 1595066ea43fSLisandro Dalcin 15969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(buf, 6, MPI_INT, 0, comm)); 1597066ea43fSLisandro Dalcin 1598066ea43fSLisandro Dalcin dim = buf[0]; 1599066ea43fSLisandro Dalcin order = buf[1]; 1600066ea43fSLisandro Dalcin periodic = buf[2] ? PETSC_TRUE : PETSC_FALSE; 1601066ea43fSLisandro Dalcin isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE; 1602066ea43fSLisandro Dalcin isHybrid = buf[4] ? PETSC_TRUE : PETSC_FALSE; 1603066ea43fSLisandro Dalcin hasTetra = buf[5] ? PETSC_TRUE : PETSC_FALSE; 1604dea2a3c7SStefano Zampini } 1605dea2a3c7SStefano Zampini 1606066ea43fSLisandro Dalcin if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE; 160708401ef6SPierre Jolivet PetscCheck(!highOrder || !isHybrid,comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet"); 1608066ea43fSLisandro Dalcin 16090598e1a2SLisandro Dalcin /* We do not want this label automatically computed, instead we fill it here */ 16109566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "celltype")); 161111c56cb0SLisandro Dalcin 1612a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 16139566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, numCells+numVerts)); 16140598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 16150598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 16160598e1a2SLisandro Dalcin DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType); 16170598e1a2SLisandro Dalcin if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR; 16189566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(*dm, cell, elem->numVerts)); 16199566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, cell, ctype)); 1620db397164SMichael Lange } 16210598e1a2SLisandro Dalcin for (v = numCells; v < numCells+numVerts; ++v) { 16229566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 162396ca5757SLisandro Dalcin } 16249566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 162596ca5757SLisandro Dalcin 1626a4bb7517SMichael Lange /* Add cell-vertex connections */ 16270598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 16280598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 16290598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 16300598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 16310598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 16320598e1a2SLisandro Dalcin cone[v] = numCells + vv; 1633db397164SMichael Lange } 16349566063dSJacob Faibussowitsch PetscCall(DMPlexReorderCell(*dm, cell, cone)); 16359566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, cell, cone)); 1636a4bb7517SMichael Lange } 163796ca5757SLisandro Dalcin 16389566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 16399566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 16409566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1641331830f3SMatthew G. Knepley if (interpolate) { 16425fd9971aSMatthew G. Knepley DM idm; 1643331830f3SMatthew G. Knepley 16449566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 16459566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1646331830f3SMatthew G. Knepley *dm = idm; 1647331830f3SMatthew G. Knepley } 1648917266a4SLisandro Dalcin 164981a1af93SMatthew G. Knepley /* Create the label "marker" over the whole boundary */ 16501dca8a05SBarry Smith PetscCheck(!usemarker || interpolate || dim <= 1,comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1651dd400576SPatrick Sanan if (rank == 0 && usemarker) { 1652d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1653d3f73514SStefano Zampini 16549566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "marker")); 16559566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd)); 1656d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1657d3f73514SStefano Zampini PetscInt suppSize; 1658d3f73514SStefano Zampini 16599566063dSJacob Faibussowitsch PetscCall(DMPlexGetSupportSize(*dm, f, &suppSize)); 1660d3f73514SStefano Zampini if (suppSize == 1) { 1661d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1662d3f73514SStefano Zampini 16639566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone)); 1664d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 16659566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &marker, "marker", cone[p], 1)); 1666d3f73514SStefano Zampini } 16679566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone)); 1668d3f73514SStefano Zampini } 1669d3f73514SStefano Zampini } 1670d3f73514SStefano Zampini } 167116ad7e67SMichael Lange 1672dd400576SPatrick Sanan if (rank == 0) { 1673a45dabc8SMatthew G. Knepley const PetscInt Nr = useregions ? mesh->numRegions : 0; 167411c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1675d1a54cefSMatthew G. Knepley 16769566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(Nr, ®ionSets)); 16779566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd)); 16780598e1a2SLisandro Dalcin for (cell = 0, e = 0; e < numElems; ++e) { 16790598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 16806e1dd89aSLawrence Mitchell 16816e1dd89aSLawrence Mitchell /* Create cell sets */ 16820598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) { 16830598e1a2SLisandro Dalcin if (elem->numTags > 0) { 1684a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1685a45dabc8SMatthew G. Knepley PetscInt r; 1686a45dabc8SMatthew G. Knepley 16879566063dSJacob Faibussowitsch if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", cell, tag)); 1688a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 16899566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], cell, tag)); 1690a45dabc8SMatthew G. Knepley } 16916e1dd89aSLawrence Mitchell } 1692917266a4SLisandro Dalcin cell++; 16936e1dd89aSLawrence Mitchell } 16940c070f12SSander Arens 16950598e1a2SLisandro Dalcin /* Create face sets */ 16960598e1a2SLisandro Dalcin if (interpolate && elem->dim == dim-1) { 16970598e1a2SLisandro Dalcin PetscInt joinSize; 16980598e1a2SLisandro Dalcin const PetscInt *join = NULL; 16995ad74b4fSMatthew G. Knepley PetscInt Nt = elem->numTags, t, r; 1700a45dabc8SMatthew G. Knepley 17010598e1a2SLisandro Dalcin /* Find the relevant facet with vertex joins */ 17020598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 17030598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 17040598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 17050598e1a2SLisandro Dalcin cone[v] = vStart + vv; 17060598e1a2SLisandro Dalcin } 17079566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join)); 170863a3b9bcSJacob 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); 17095ad74b4fSMatthew G. Knepley for (t = 0; t < Nt; ++t) { 17105ad74b4fSMatthew G. Knepley const PetscInt tag = elem->tags[t]; 17115ad74b4fSMatthew G. Knepley 17129566063dSJacob Faibussowitsch if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], tag)); 1713a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 17149566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], join[0], tag)); 1715a45dabc8SMatthew G. Knepley } 17165ad74b4fSMatthew G. Knepley } 17179566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join)); 17180598e1a2SLisandro Dalcin } 17190598e1a2SLisandro Dalcin 17200c070f12SSander Arens /* Create vertex sets */ 17210598e1a2SLisandro Dalcin if (elem->dim == 0) { 17220598e1a2SLisandro Dalcin if (elem->numTags > 0) { 17230598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[0]; 17240598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1725a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1726a45dabc8SMatthew G. Knepley PetscInt r; 1727a45dabc8SMatthew G. Knepley 17289566063dSJacob Faibussowitsch if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag)); 172981a1af93SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 17309566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag)); 173181a1af93SMatthew G. Knepley } 173281a1af93SMatthew G. Knepley } 173381a1af93SMatthew G. Knepley } 173481a1af93SMatthew G. Knepley } 173581a1af93SMatthew G. Knepley if (markvertices) { 173681a1af93SMatthew G. Knepley for (v = 0; v < numNodes; ++v) { 173781a1af93SMatthew G. Knepley const PetscInt vv = mesh->vertexMap[v]; 17387d5b9244SMatthew G. Knepley const PetscInt *tags = &mesh->nodelist->tag[v*GMSH_MAX_TAGS]; 17397d5b9244SMatthew G. Knepley PetscInt r, t; 174081a1af93SMatthew G. Knepley 17417d5b9244SMatthew G. Knepley for (t = 0; t < GMSH_MAX_TAGS; ++t) { 17427d5b9244SMatthew G. Knepley const PetscInt tag = tags[t]; 17437d5b9244SMatthew G. Knepley 17447d5b9244SMatthew G. Knepley if (tag == -1) continue; 17459566063dSJacob Faibussowitsch if (!Nr) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag)); 1746a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 17479566063dSJacob Faibussowitsch if (mesh->regionTags[r] == tag) PetscCall(DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag)); 17480598e1a2SLisandro Dalcin } 17490598e1a2SLisandro Dalcin } 17500598e1a2SLisandro Dalcin } 17510598e1a2SLisandro Dalcin } 17529566063dSJacob Faibussowitsch PetscCall(PetscFree(regionSets)); 1753a45dabc8SMatthew G. Knepley } 17540598e1a2SLisandro Dalcin 17557dd454faSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 17567dd454faSLisandro Dalcin enum {n = 4}; 17577dd454faSLisandro Dalcin PetscBool flag[n]; 17587dd454faSLisandro Dalcin 17597dd454faSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 17607dd454faSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 17617dd454faSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 17627dd454faSLisandro Dalcin flag[3] = marker ? PETSC_TRUE : PETSC_FALSE; 17639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 17649566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 17659566063dSJacob Faibussowitsch if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 17669566063dSJacob Faibussowitsch if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 17679566063dSJacob Faibussowitsch if (flag[3]) PetscCall(DMCreateLabel(*dm, "marker")); 17687dd454faSLisandro Dalcin } 17697dd454faSLisandro Dalcin 17700598e1a2SLisandro Dalcin if (periodic) { 17719566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(numVerts, &periodicVerts)); 17720598e1a2SLisandro Dalcin for (n = 0; n < numNodes; ++n) { 17730598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) { 17740598e1a2SLisandro Dalcin if (PetscUnlikely(mesh->periodMap[n] != n)) { 17750598e1a2SLisandro Dalcin PetscInt m = mesh->periodMap[n]; 17769566063dSJacob Faibussowitsch PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[n])); 17779566063dSJacob Faibussowitsch PetscCall(PetscBTSet(periodicVerts, mesh->vertexMap[m])); 17780598e1a2SLisandro Dalcin } 17790598e1a2SLisandro Dalcin } 17800598e1a2SLisandro Dalcin } 17819566063dSJacob Faibussowitsch PetscCall(PetscBTCreate(numCells, &periodicCells)); 17820598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 17830598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 17840598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 17850598e1a2SLisandro Dalcin PetscInt nn = elem->nodes[v]; 17860598e1a2SLisandro Dalcin PetscInt vv = mesh->vertexMap[nn]; 17870598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) { 17889566063dSJacob Faibussowitsch PetscCall(PetscBTSet(periodicCells, cell)); break; 17890c070f12SSander Arens } 17900c070f12SSander Arens } 17910c070f12SSander Arens } 179216ad7e67SMichael Lange } 179316ad7e67SMichael Lange 1794066ea43fSLisandro Dalcin /* Setup coordinate DM */ 17950598e1a2SLisandro Dalcin if (coordDim < 0) coordDim = dim; 17969566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, coordDim)); 17979566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(*dm, &cdm)); 1798066ea43fSLisandro Dalcin if (highOrder) { 1799066ea43fSLisandro Dalcin PetscFE fe; 1800066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1801066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_EQUISPACED; 1802066ea43fSLisandro Dalcin 1803066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1804066ea43fSLisandro Dalcin 18059566063dSJacob Faibussowitsch PetscCall(GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe)); 18069566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view")); 18079566063dSJacob Faibussowitsch PetscCall(DMSetField(cdm, 0, NULL, (PetscObject) fe)); 18089566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 18099566063dSJacob Faibussowitsch PetscCall(DMCreateDS(cdm)); 1810066ea43fSLisandro Dalcin } 18116858538eSMatthew G. Knepley if (periodic) { 18126858538eSMatthew G. Knepley PetscCall(DMClone(cdm, &cdmCell)); 18136858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateDM(*dm, cdmCell)); 18146858538eSMatthew G. Knepley } 1815066ea43fSLisandro Dalcin 1816066ea43fSLisandro Dalcin /* Create coordinates */ 1817066ea43fSLisandro Dalcin if (highOrder) { 1818066ea43fSLisandro Dalcin PetscInt maxDof = GmshNumNodes_HEX(order)*coordDim; 1819066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1820066ea43fSLisandro Dalcin PetscSection section; 1821066ea43fSLisandro Dalcin PetscScalar *cellCoords; 1822066ea43fSLisandro Dalcin 18239566063dSJacob Faibussowitsch PetscCall(DMSetLocalSection(cdm, NULL)); 18246858538eSMatthew G. Knepley PetscCall(DMGetLocalSection(cdm, &cs)); 18256858538eSMatthew G. Knepley PetscCall(PetscSectionClone(cs, §ion)); 18269566063dSJacob Faibussowitsch PetscCall(DMPlexSetClosurePermutationTensor(cdm, 0, section)); /* XXX Implement DMPlexSetClosurePermutationLexicographic() */ 1827066ea43fSLisandro Dalcin 18289566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(cdm, &coordinates)); 18299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(maxDof, &cellCoords)); 1830066ea43fSLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1831066ea43fSLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1832066ea43fSLisandro Dalcin const int *lexorder = GmshCellMap[elem->cellType].lexorder(); 1833b9bf55e5SMatthew G. Knepley int s = 0; 1834066ea43fSLisandro Dalcin for (n = 0; n < elem->numNodes; ++n) { 1835b9bf55e5SMatthew G. Knepley while (lexorder[n+s] < 0) ++s; 1836b9bf55e5SMatthew G. Knepley const PetscInt node = elem->nodes[lexorder[n+s]]; 1837b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[(n+s)*coordDim+d] = (PetscReal) coords[node*3+d]; 1838b9bf55e5SMatthew G. Knepley } 1839b9bf55e5SMatthew G. Knepley if (s) { 1840b9bf55e5SMatthew G. Knepley /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */ 1841b9bf55e5SMatthew G. Knepley PetscReal quaCenterWeights[9] = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25}; 1842b9bf55e5SMatthew G. Knepley /* For the coordinate mapping we weight vertices by -1/4 and edges by 1/2, which we get from Q_2 interpolation */ 1843b9bf55e5SMatthew G. Knepley PetscReal hexBottomWeights[27] = {-0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25, 1844b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1845b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}; 1846b9bf55e5SMatthew G. Knepley PetscReal hexFrontWeights[27] = {-0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1847b9bf55e5SMatthew G. Knepley 0.5, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1848b9bf55e5SMatthew G. Knepley -0.25, 0.5, -0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0}; 1849b9bf55e5SMatthew G. Knepley PetscReal hexLeftWeights[27] = {-0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0, 1850b9bf55e5SMatthew G. Knepley 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 1851b9bf55e5SMatthew G. Knepley -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 0.0, 0.0}; 1852b9bf55e5SMatthew G. Knepley PetscReal hexRightWeights[27] = { 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25, 1853b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1854b9bf55e5SMatthew G. Knepley 0.0, 0.0, -0.25, 0.0, 0.0, 0.5, 0.0, 0.0, -0.25}; 1855b9bf55e5SMatthew G. Knepley PetscReal hexBackWeights[27] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25, 1856b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.5, 1857b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, -0.25, 0.5, -0.25}; 1858b9bf55e5SMatthew G. Knepley PetscReal hexTopWeights[27] = { 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1859b9bf55e5SMatthew G. Knepley 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1860b9bf55e5SMatthew G. Knepley -0.25, 0.5, -0.25, 0.5, 0.0, 0.5, -0.25, 0.5, -0.25}; 1861b9bf55e5SMatthew G. Knepley PetscReal hexCenterWeights[27] = {-0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25, 1862b9bf55e5SMatthew G. Knepley 0.25, 0.0, 0.25, 0.0, 0.0, 0.0, 0.25, 0.0, 0.25, 1863b9bf55e5SMatthew G. Knepley -0.25, 0.25, -0.25, 0.25, 0.0, 0.25, -0.25, 0.25, -0.25}; 1864b9bf55e5SMatthew G. Knepley PetscReal *sdWeights2[9] = {NULL, NULL, NULL, NULL, quaCenterWeights, NULL, NULL, NULL, NULL}; 1865b9bf55e5SMatthew G. Knepley PetscReal *sdWeights3[27] = {NULL, NULL, NULL, NULL, hexBottomWeights, NULL, NULL, NULL, NULL, 1866b9bf55e5SMatthew G. Knepley NULL, hexFrontWeights, NULL, hexLeftWeights, hexCenterWeights, hexRightWeights, NULL, hexBackWeights, NULL, 1867b9bf55e5SMatthew G. Knepley NULL, NULL, NULL, NULL, hexTopWeights, NULL, NULL, NULL, NULL}; 1868b9bf55e5SMatthew G. Knepley PetscReal **sdWeights[4] = {NULL, NULL, sdWeights2, sdWeights3}; 1869b9bf55e5SMatthew G. Knepley 1870b9bf55e5SMatthew G. Knepley /* Missing entries in serendipity cell, only works for 8-node quad and 20-node hex */ 1871b9bf55e5SMatthew G. Knepley for (n = 0; n < elem->numNodes+s; ++n) { 1872b9bf55e5SMatthew G. Knepley if (lexorder[n] >= 0) continue; 1873b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[n*coordDim+d] = 0.0; 1874b9bf55e5SMatthew G. Knepley for (int bn = 0; bn < elem->numNodes+s; ++bn) { 1875b9bf55e5SMatthew G. Knepley if (lexorder[bn] < 0) continue; 1876b9bf55e5SMatthew G. Knepley const PetscReal *weights = sdWeights[coordDim][n]; 1877b9bf55e5SMatthew G. Knepley const PetscInt bnode = elem->nodes[lexorder[bn]]; 1878b9bf55e5SMatthew G. Knepley for (d = 0; d < coordDim; ++d) cellCoords[n*coordDim+d] += weights[bn] * (PetscReal) coords[bnode*3+d]; 1879b9bf55e5SMatthew G. Knepley } 1880b9bf55e5SMatthew G. Knepley } 1881066ea43fSLisandro Dalcin } 18829566063dSJacob Faibussowitsch PetscCall(DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES)); 1883066ea43fSLisandro Dalcin } 18849566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(§ion)); 18859566063dSJacob Faibussowitsch PetscCall(PetscFree(cellCoords)); 1886066ea43fSLisandro Dalcin } else { 1887066ea43fSLisandro Dalcin PetscInt *nodeMap; 1888066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1889066ea43fSLisandro Dalcin PetscScalar *pointCoords; 1890066ea43fSLisandro Dalcin 18916858538eSMatthew G. Knepley PetscCall(DMGetCoordinateSection(*dm, &cs)); 18926858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(cs, 1)); 18936858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(cs, 0, coordDim)); 18946858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(cs, numCells, numCells+numVerts)); 18956858538eSMatthew G. Knepley for (v = numCells; v < numCells+numVerts; ++v) { 18966858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(cs, v, coordDim)); 18976858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(cs, v, 0, coordDim)); 1898f45c9589SStefano Zampini } 18996858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(cs)); 19006858538eSMatthew G. Knepley 19016858538eSMatthew G. Knepley /* We need to localize coordinates on cells */ 19020598e1a2SLisandro Dalcin if (periodic) { 19036858538eSMatthew G. Knepley PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) cdmCell), &csCell)); 19046858538eSMatthew G. Knepley PetscCall(PetscSectionSetNumFields(csCell, 1)); 19056858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldComponents(csCell, 0, coordDim)); 19066858538eSMatthew G. Knepley PetscCall(PetscSectionSetChart(csCell, 0, numCells)); 19070598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 19080598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 19090598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 19100598e1a2SLisandro Dalcin PetscInt dof = elem->numVerts * coordDim; 19116858538eSMatthew G. Knepley 19126858538eSMatthew G. Knepley PetscCall(PetscSectionSetDof(csCell, cell, dof)); 19136858538eSMatthew G. Knepley PetscCall(PetscSectionSetFieldDof(csCell, cell, 0, dof)); 1914f45c9589SStefano Zampini } 1915f45c9589SStefano Zampini } 19166858538eSMatthew G. Knepley PetscCall(PetscSectionSetUp(csCell)); 19176858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinateSection(*dm, PETSC_DETERMINE, csCell)); 1918f45c9589SStefano Zampini } 1919066ea43fSLisandro Dalcin 19209566063dSJacob Faibussowitsch PetscCall(DMCreateLocalVector(cdm, &coordinates)); 19219566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &pointCoords)); 19226858538eSMatthew G. Knepley PetscCall(PetscMalloc1(numVerts, &nodeMap)); 19236858538eSMatthew G. Knepley for (n = 0; n < numNodes; n++) 19246858538eSMatthew G. Knepley if (mesh->vertexMap[n] >= 0) nodeMap[mesh->vertexMap[n]] = n; 19256858538eSMatthew G. Knepley for (v = 0; v < numVerts; ++v) { 19266858538eSMatthew G. Knepley PetscInt off, node = nodeMap[v]; 19276858538eSMatthew G. Knepley 19286858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(cs, numCells + v, &off)); 19296858538eSMatthew G. Knepley for (d = 0; d < coordDim; ++d) pointCoords[off+d] = (PetscReal) coords[node*3+d]; 19306858538eSMatthew G. Knepley } 19316858538eSMatthew G. Knepley PetscCall(VecRestoreArray(coordinates, &pointCoords)); 19326858538eSMatthew G. Knepley PetscCall(PetscFree(nodeMap)); 19336858538eSMatthew G. Knepley 19340598e1a2SLisandro Dalcin if (periodic) { 19356858538eSMatthew G. Knepley PetscCall(DMCreateLocalVector(cdmCell, &coordinatesCell)); 19366858538eSMatthew G. Knepley PetscCall(VecGetArray(coordinatesCell, &pointCoords)); 19370598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 19380598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 19390598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 19400598e1a2SLisandro Dalcin PetscInt off, node; 19416858538eSMatthew G. Knepley for (v = 0; v < elem->numVerts; ++v) cone[v] = elem->nodes[v]; 19426858538eSMatthew G. Knepley PetscCall(DMPlexReorderCell(cdmCell, cell, cone)); 19436858538eSMatthew G. Knepley PetscCall(PetscSectionGetOffset(csCell, cell, &off)); 19440598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 19450598e1a2SLisandro Dalcin for (node = cone[v], d = 0; d < coordDim; ++d) 1946066ea43fSLisandro Dalcin pointCoords[off++] = (PetscReal) coords[node*3+d]; 1947331830f3SMatthew G. Knepley } 1948331830f3SMatthew G. Knepley } 19496858538eSMatthew G. Knepley PetscCall(VecSetBlockSize(coordinatesCell, coordDim)); 19506858538eSMatthew G. Knepley PetscCall(VecRestoreArray(coordinatesCell, &pointCoords)); 19516858538eSMatthew G. Knepley PetscCall(DMSetCellCoordinatesLocal(*dm, coordinatesCell)); 19526858538eSMatthew G. Knepley PetscCall(VecDestroy(&coordinatesCell)); 1953331830f3SMatthew G. Knepley } 19546858538eSMatthew G. Knepley PetscCall(PetscSectionDestroy(&csCell)); 19556858538eSMatthew G. Knepley PetscCall(DMDestroy(&cdmCell)); 1956066ea43fSLisandro Dalcin } 1957066ea43fSLisandro Dalcin 19589566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 19599566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, coordDim)); 19609566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 19619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 196211c56cb0SLisandro Dalcin 19639566063dSJacob Faibussowitsch PetscCall(GmshMeshDestroy(&mesh)); 19649566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&periodicVerts)); 19659566063dSJacob Faibussowitsch PetscCall(PetscBTDestroy(&periodicCells)); 196611c56cb0SLisandro Dalcin 1967066ea43fSLisandro Dalcin if (highOrder && project) { 1968066ea43fSLisandro Dalcin PetscFE fe; 1969066ea43fSLisandro Dalcin const char prefix[] = "dm_plex_gmsh_project_"; 1970066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1971066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI; 1972066ea43fSLisandro Dalcin 1973066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1974066ea43fSLisandro Dalcin 19759566063dSJacob Faibussowitsch PetscCall(GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe)); 19769566063dSJacob Faibussowitsch PetscCall(PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view")); 19779566063dSJacob Faibussowitsch PetscCall(DMProjectCoordinates(*dm, fe)); 19789566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 1979066ea43fSLisandro Dalcin } 1980066ea43fSLisandro Dalcin 19819566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL)); 1982331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1983331830f3SMatthew G. Knepley } 1984