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 16066ea43fSLisandro Dalcin #define GMSH_LEXORDER_LIST(T) \ 17066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 1) \ 18066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 2) \ 19066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 3) \ 20066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 4) \ 21066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 5) \ 22066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 6) \ 23066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 7) \ 24066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 8) \ 25066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 9) \ 26066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(T, 10) 27066ea43fSLisandro Dalcin 28066ea43fSLisandro Dalcin GMSH_LEXORDER_ITEM(VTX, 0) 29066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(SEG) 30066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TRI) 31066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(QUA) 32066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(TET) 33066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(HEX) 34066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PRI) 35066ea43fSLisandro Dalcin GMSH_LEXORDER_LIST(PYR) 36066ea43fSLisandro Dalcin 37066ea43fSLisandro Dalcin typedef enum { 38066ea43fSLisandro Dalcin GMSH_VTX = 0, 39066ea43fSLisandro Dalcin GMSH_SEG = 1, 40066ea43fSLisandro Dalcin GMSH_TRI = 2, 41066ea43fSLisandro Dalcin GMSH_QUA = 3, 42066ea43fSLisandro Dalcin GMSH_TET = 4, 43066ea43fSLisandro Dalcin GMSH_HEX = 5, 44066ea43fSLisandro Dalcin GMSH_PRI = 6, 45066ea43fSLisandro Dalcin GMSH_PYR = 7, 46066ea43fSLisandro Dalcin GMSH_NUM_POLYTOPES = 8 47066ea43fSLisandro Dalcin } GmshPolytopeType; 48066ea43fSLisandro Dalcin 49de78e4feSLisandro Dalcin typedef struct { 500598e1a2SLisandro Dalcin int cellType; 51066ea43fSLisandro Dalcin int polytope; 520598e1a2SLisandro Dalcin int dim; 530598e1a2SLisandro Dalcin int order; 54066ea43fSLisandro Dalcin int numVerts; 550598e1a2SLisandro Dalcin int numNodes; 56066ea43fSLisandro Dalcin int* (*lexorder)(void); 570598e1a2SLisandro Dalcin } GmshCellInfo; 580598e1a2SLisandro Dalcin 59066ea43fSLisandro Dalcin #define GmshCellEntry(cellType, polytope, dim, order) \ 60066ea43fSLisandro Dalcin {cellType, GMSH_##polytope, dim, order, \ 61066ea43fSLisandro Dalcin GmshNumNodes_##polytope(1), \ 62066ea43fSLisandro Dalcin GmshNumNodes_##polytope(order), \ 63066ea43fSLisandro Dalcin GmshLexOrder_##polytope##_##order} 640598e1a2SLisandro Dalcin 650598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = { 66066ea43fSLisandro Dalcin GmshCellEntry( 15, VTX, 0, 0), 670598e1a2SLisandro Dalcin 68066ea43fSLisandro Dalcin GmshCellEntry( 1, SEG, 1, 1), 69066ea43fSLisandro Dalcin GmshCellEntry( 8, SEG, 1, 2), 70066ea43fSLisandro Dalcin GmshCellEntry( 26, SEG, 1, 3), 71066ea43fSLisandro Dalcin GmshCellEntry( 27, SEG, 1, 4), 72066ea43fSLisandro Dalcin GmshCellEntry( 28, SEG, 1, 5), 73066ea43fSLisandro Dalcin GmshCellEntry( 62, SEG, 1, 6), 74066ea43fSLisandro Dalcin GmshCellEntry( 63, SEG, 1, 7), 75066ea43fSLisandro Dalcin GmshCellEntry( 64, SEG, 1, 8), 76066ea43fSLisandro Dalcin GmshCellEntry( 65, SEG, 1, 9), 77066ea43fSLisandro Dalcin GmshCellEntry( 66, SEG, 1, 10), 780598e1a2SLisandro Dalcin 79066ea43fSLisandro Dalcin GmshCellEntry( 2, TRI, 2, 1), 80066ea43fSLisandro Dalcin GmshCellEntry( 9, TRI, 2, 2), 81066ea43fSLisandro Dalcin GmshCellEntry( 21, TRI, 2, 3), 82066ea43fSLisandro Dalcin GmshCellEntry( 23, TRI, 2, 4), 83066ea43fSLisandro Dalcin GmshCellEntry( 25, TRI, 2, 5), 84066ea43fSLisandro Dalcin GmshCellEntry( 42, TRI, 2, 6), 85066ea43fSLisandro Dalcin GmshCellEntry( 43, TRI, 2, 7), 86066ea43fSLisandro Dalcin GmshCellEntry( 44, TRI, 2, 8), 87066ea43fSLisandro Dalcin GmshCellEntry( 45, TRI, 2, 9), 88066ea43fSLisandro Dalcin GmshCellEntry( 46, TRI, 2, 10), 890598e1a2SLisandro Dalcin 90066ea43fSLisandro Dalcin GmshCellEntry( 3, QUA, 2, 1), 91066ea43fSLisandro Dalcin GmshCellEntry( 10, QUA, 2, 2), 92066ea43fSLisandro Dalcin GmshCellEntry( 36, QUA, 2, 3), 93066ea43fSLisandro Dalcin GmshCellEntry( 37, QUA, 2, 4), 94066ea43fSLisandro Dalcin GmshCellEntry( 38, QUA, 2, 5), 95066ea43fSLisandro Dalcin GmshCellEntry( 47, QUA, 2, 6), 96066ea43fSLisandro Dalcin GmshCellEntry( 48, QUA, 2, 7), 97066ea43fSLisandro Dalcin GmshCellEntry( 49, QUA, 2, 8), 98066ea43fSLisandro Dalcin GmshCellEntry( 50, QUA, 2, 9), 99066ea43fSLisandro Dalcin GmshCellEntry( 51, QUA, 2, 10), 1000598e1a2SLisandro Dalcin 101066ea43fSLisandro Dalcin GmshCellEntry( 4, TET, 3, 1), 102066ea43fSLisandro Dalcin GmshCellEntry( 11, TET, 3, 2), 103066ea43fSLisandro Dalcin GmshCellEntry( 29, TET, 3, 3), 104066ea43fSLisandro Dalcin GmshCellEntry( 30, TET, 3, 4), 105066ea43fSLisandro Dalcin GmshCellEntry( 31, TET, 3, 5), 106066ea43fSLisandro Dalcin GmshCellEntry( 71, TET, 3, 6), 107066ea43fSLisandro Dalcin GmshCellEntry( 72, TET, 3, 7), 108066ea43fSLisandro Dalcin GmshCellEntry( 73, TET, 3, 8), 109066ea43fSLisandro Dalcin GmshCellEntry( 74, TET, 3, 9), 110066ea43fSLisandro Dalcin GmshCellEntry( 75, TET, 3, 10), 1110598e1a2SLisandro Dalcin 112066ea43fSLisandro Dalcin GmshCellEntry( 5, HEX, 3, 1), 113066ea43fSLisandro Dalcin GmshCellEntry( 12, HEX, 3, 2), 114066ea43fSLisandro Dalcin GmshCellEntry( 92, HEX, 3, 3), 115066ea43fSLisandro Dalcin GmshCellEntry( 93, HEX, 3, 4), 116066ea43fSLisandro Dalcin GmshCellEntry( 94, HEX, 3, 5), 117066ea43fSLisandro Dalcin GmshCellEntry( 95, HEX, 3, 6), 118066ea43fSLisandro Dalcin GmshCellEntry( 96, HEX, 3, 7), 119066ea43fSLisandro Dalcin GmshCellEntry( 97, HEX, 3, 8), 120066ea43fSLisandro Dalcin GmshCellEntry( 98, HEX, 3, 9), 121066ea43fSLisandro Dalcin GmshCellEntry( -1, HEX, 3, 10), 1220598e1a2SLisandro Dalcin 123066ea43fSLisandro Dalcin GmshCellEntry( 6, PRI, 3, 1), 124066ea43fSLisandro Dalcin GmshCellEntry( 13, PRI, 3, 2), 125066ea43fSLisandro Dalcin GmshCellEntry( 90, PRI, 3, 3), 126066ea43fSLisandro Dalcin GmshCellEntry( 91, PRI, 3, 4), 127066ea43fSLisandro Dalcin GmshCellEntry(106, PRI, 3, 5), 128066ea43fSLisandro Dalcin GmshCellEntry(107, PRI, 3, 6), 129066ea43fSLisandro Dalcin GmshCellEntry(108, PRI, 3, 7), 130066ea43fSLisandro Dalcin GmshCellEntry(109, PRI, 3, 8), 131066ea43fSLisandro Dalcin GmshCellEntry(110, PRI, 3, 9), 132066ea43fSLisandro Dalcin GmshCellEntry( -1, PRI, 3, 10), 1330598e1a2SLisandro Dalcin 134066ea43fSLisandro Dalcin GmshCellEntry( 7, PYR, 3, 1), 135066ea43fSLisandro Dalcin GmshCellEntry( 14, PYR, 3, 2), 136066ea43fSLisandro Dalcin GmshCellEntry(118, PYR, 3, 3), 137066ea43fSLisandro Dalcin GmshCellEntry(119, PYR, 3, 4), 138066ea43fSLisandro Dalcin GmshCellEntry(120, PYR, 3, 5), 139066ea43fSLisandro Dalcin GmshCellEntry(121, PYR, 3, 6), 140066ea43fSLisandro Dalcin GmshCellEntry(122, PYR, 3, 7), 141066ea43fSLisandro Dalcin GmshCellEntry(123, PYR, 3, 8), 142066ea43fSLisandro Dalcin GmshCellEntry(124, PYR, 3, 9), 143066ea43fSLisandro Dalcin GmshCellEntry( -1, PYR, 3, 10) 1440598e1a2SLisandro Dalcin 1450598e1a2SLisandro Dalcin #if 0 146066ea43fSLisandro Dalcin {20, GMSH_TRI, 2, 3, 3, 9, NULL}, 147066ea43fSLisandro Dalcin {16, GMSH_QUA, 2, 2, 4, 8, NULL}, 148066ea43fSLisandro Dalcin {17, GMSH_HEX, 3, 2, 8, 20, NULL}, 149066ea43fSLisandro Dalcin {18, GMSH_PRI, 3, 2, 6, 15, NULL}, 150066ea43fSLisandro Dalcin {19, GMSH_PYR, 3, 2, 5, 13, NULL}, 1510598e1a2SLisandro Dalcin #endif 1520598e1a2SLisandro Dalcin }; 1530598e1a2SLisandro Dalcin 1540598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150]; 1550598e1a2SLisandro Dalcin 1560598e1a2SLisandro Dalcin static PetscErrorCode GmshCellInfoSetUp(void) 1570598e1a2SLisandro Dalcin { 1580598e1a2SLisandro Dalcin size_t i, n; 1590598e1a2SLisandro Dalcin static PetscBool called = PETSC_FALSE; 1600598e1a2SLisandro Dalcin 1610598e1a2SLisandro Dalcin if (called) return 0; 1620598e1a2SLisandro Dalcin PetscFunctionBegin; 1630598e1a2SLisandro Dalcin called = PETSC_TRUE; 1640598e1a2SLisandro Dalcin n = sizeof(GmshCellMap)/sizeof(GmshCellMap[0]); 1650598e1a2SLisandro Dalcin for (i = 0; i < n; ++i) { 1660598e1a2SLisandro Dalcin GmshCellMap[i].cellType = -1; 167066ea43fSLisandro Dalcin GmshCellMap[i].polytope = -1; 1680598e1a2SLisandro Dalcin } 1690598e1a2SLisandro Dalcin n = sizeof(GmshCellTable)/sizeof(GmshCellTable[0]); 170066ea43fSLisandro Dalcin for (i = 0; i < n; ++i) { 171066ea43fSLisandro Dalcin if (GmshCellTable[i].cellType <= 0) continue; 172066ea43fSLisandro Dalcin GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i]; 173066ea43fSLisandro Dalcin } 1740598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1750598e1a2SLisandro Dalcin } 1760598e1a2SLisandro Dalcin 1770598e1a2SLisandro Dalcin #define GmshCellTypeCheck(ct) 0; do { \ 1780598e1a2SLisandro Dalcin const int _ct_ = (int)ct; \ 1790598e1a2SLisandro Dalcin if (_ct_ < 0 || _ct_ >= (int)(sizeof(GmshCellMap)/sizeof(GmshCellMap[0]))) \ 1800598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); \ 1810598e1a2SLisandro Dalcin if (GmshCellMap[_ct_].cellType != _ct_) \ 1820598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 183066ea43fSLisandro Dalcin if (GmshCellMap[_ct_].polytope == -1) \ 1840598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 1850598e1a2SLisandro Dalcin } while (0) 1860598e1a2SLisandro Dalcin 1870598e1a2SLisandro Dalcin typedef struct { 188698a79b9SLisandro Dalcin PetscViewer viewer; 189698a79b9SLisandro Dalcin int fileFormat; 190698a79b9SLisandro Dalcin int dataSize; 191698a79b9SLisandro Dalcin PetscBool binary; 192698a79b9SLisandro Dalcin PetscBool byteSwap; 193698a79b9SLisandro Dalcin size_t wlen; 194698a79b9SLisandro Dalcin void *wbuf; 195698a79b9SLisandro Dalcin size_t slen; 196698a79b9SLisandro Dalcin void *sbuf; 1970598e1a2SLisandro Dalcin PetscInt *nbuf; 1980598e1a2SLisandro Dalcin PetscInt nodeStart; 1990598e1a2SLisandro Dalcin PetscInt nodeEnd; 20033a088b6SMatthew G. Knepley PetscInt *nodeMap; 201698a79b9SLisandro Dalcin } GmshFile; 202de78e4feSLisandro Dalcin 203698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) 204de78e4feSLisandro Dalcin { 205de78e4feSLisandro Dalcin size_t size = count * eltsize; 20611c56cb0SLisandro Dalcin PetscErrorCode ierr; 20711c56cb0SLisandro Dalcin 20811c56cb0SLisandro Dalcin PetscFunctionBegin; 209698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 210698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 211698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr); 212698a79b9SLisandro Dalcin gmsh->wlen = size; 213de78e4feSLisandro Dalcin } 214698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->wbuf : NULL; 215de78e4feSLisandro Dalcin PetscFunctionReturn(0); 216de78e4feSLisandro Dalcin } 217de78e4feSLisandro Dalcin 218698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) 219698a79b9SLisandro Dalcin { 220698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 221698a79b9SLisandro Dalcin size_t size = count * dataSize; 222698a79b9SLisandro Dalcin PetscErrorCode ierr; 223698a79b9SLisandro Dalcin 224698a79b9SLisandro Dalcin PetscFunctionBegin; 225698a79b9SLisandro Dalcin if (gmsh->slen < size) { 226698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 227698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr); 228698a79b9SLisandro Dalcin gmsh->slen = size; 229698a79b9SLisandro Dalcin } 230698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->sbuf : NULL; 231698a79b9SLisandro Dalcin PetscFunctionReturn(0); 232698a79b9SLisandro Dalcin } 233698a79b9SLisandro Dalcin 234698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) 235de78e4feSLisandro Dalcin { 236de78e4feSLisandro Dalcin PetscErrorCode ierr; 237de78e4feSLisandro Dalcin PetscFunctionBegin; 238698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr); 239698a79b9SLisandro Dalcin if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);} 240698a79b9SLisandro Dalcin PetscFunctionReturn(0); 241698a79b9SLisandro Dalcin } 242698a79b9SLisandro Dalcin 243698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) 244698a79b9SLisandro Dalcin { 245698a79b9SLisandro Dalcin PetscErrorCode ierr; 246698a79b9SLisandro Dalcin PetscFunctionBegin; 247698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr); 248698a79b9SLisandro Dalcin PetscFunctionReturn(0); 249698a79b9SLisandro Dalcin } 250698a79b9SLisandro Dalcin 251d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) 252d6689ca9SLisandro Dalcin { 253d6689ca9SLisandro Dalcin PetscErrorCode ierr; 254d6689ca9SLisandro Dalcin PetscFunctionBegin; 255d6689ca9SLisandro Dalcin ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr); 256d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 257d6689ca9SLisandro Dalcin } 258d6689ca9SLisandro Dalcin 259d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) 260d6689ca9SLisandro Dalcin { 261d6689ca9SLisandro Dalcin PetscBool match; 262d6689ca9SLisandro Dalcin PetscErrorCode ierr; 263d6689ca9SLisandro Dalcin 264d6689ca9SLisandro Dalcin PetscFunctionBegin; 265d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr); 2660598e1a2SLisandro Dalcin if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s",Section); 267d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 268d6689ca9SLisandro Dalcin } 269d6689ca9SLisandro Dalcin 270d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) 271d6689ca9SLisandro Dalcin { 272d6689ca9SLisandro Dalcin PetscBool match; 273d6689ca9SLisandro Dalcin PetscErrorCode ierr; 274d6689ca9SLisandro Dalcin 275d6689ca9SLisandro Dalcin PetscFunctionBegin; 276d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 277d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 278d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr); 279d6689ca9SLisandro Dalcin if (!match) break; 280d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 281d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 282d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr); 283d6689ca9SLisandro Dalcin if (match) break; 284d6689ca9SLisandro Dalcin } 285d6689ca9SLisandro Dalcin } 286d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 287d6689ca9SLisandro Dalcin } 288d6689ca9SLisandro Dalcin 289d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) 290d6689ca9SLisandro Dalcin { 291d6689ca9SLisandro Dalcin PetscErrorCode ierr; 292d6689ca9SLisandro Dalcin PetscFunctionBegin; 293d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 294d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr); 295d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 296d6689ca9SLisandro Dalcin } 297d6689ca9SLisandro Dalcin 298698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) 299698a79b9SLisandro Dalcin { 300698a79b9SLisandro Dalcin PetscInt i; 301698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 302698a79b9SLisandro Dalcin PetscErrorCode ierr; 303698a79b9SLisandro Dalcin 304698a79b9SLisandro Dalcin PetscFunctionBegin; 305698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 306698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr); 307698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 308698a79b9SLisandro Dalcin int *ibuf = NULL; 30989d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 310698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 311698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 312698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 313698a79b9SLisandro Dalcin long *ibuf = NULL; 31489d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 315698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr); 316698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 317698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 318698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 31989d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 320698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr); 321698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 322698a79b9SLisandro Dalcin } 323698a79b9SLisandro Dalcin PetscFunctionReturn(0); 324698a79b9SLisandro Dalcin } 325698a79b9SLisandro Dalcin 326698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) 327698a79b9SLisandro Dalcin { 328698a79b9SLisandro Dalcin PetscErrorCode ierr; 329698a79b9SLisandro Dalcin PetscFunctionBegin; 330698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr); 331698a79b9SLisandro Dalcin PetscFunctionReturn(0); 332698a79b9SLisandro Dalcin } 333698a79b9SLisandro Dalcin 334698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) 335698a79b9SLisandro Dalcin { 336698a79b9SLisandro Dalcin PetscErrorCode ierr; 337698a79b9SLisandro Dalcin PetscFunctionBegin; 338698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr); 339de78e4feSLisandro Dalcin PetscFunctionReturn(0); 340de78e4feSLisandro Dalcin } 341de78e4feSLisandro Dalcin 342de78e4feSLisandro Dalcin typedef struct { 3430598e1a2SLisandro Dalcin PetscInt id; /* Entity ID */ 3440598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 345de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 346de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 347de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 348de78e4feSLisandro Dalcin } GmshEntity; 349de78e4feSLisandro Dalcin 350de78e4feSLisandro Dalcin typedef struct { 351730356f6SLisandro Dalcin GmshEntity *entity[4]; 352730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 353730356f6SLisandro Dalcin } GmshEntities; 354730356f6SLisandro Dalcin 355698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) 356730356f6SLisandro Dalcin { 357730356f6SLisandro Dalcin PetscInt dim; 358730356f6SLisandro Dalcin PetscErrorCode ierr; 359730356f6SLisandro Dalcin 360730356f6SLisandro Dalcin PetscFunctionBegin; 361730356f6SLisandro Dalcin ierr = PetscNew(entities);CHKERRQ(ierr); 362730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 363730356f6SLisandro Dalcin ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr); 364730356f6SLisandro Dalcin ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 365730356f6SLisandro Dalcin } 366730356f6SLisandro Dalcin PetscFunctionReturn(0); 367730356f6SLisandro Dalcin } 368730356f6SLisandro Dalcin 3690598e1a2SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 3700598e1a2SLisandro Dalcin { 3710598e1a2SLisandro Dalcin PetscInt dim; 3720598e1a2SLisandro Dalcin PetscErrorCode ierr; 3730598e1a2SLisandro Dalcin 3740598e1a2SLisandro Dalcin PetscFunctionBegin; 3750598e1a2SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 3760598e1a2SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 3770598e1a2SLisandro Dalcin ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr); 3780598e1a2SLisandro Dalcin ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 3790598e1a2SLisandro Dalcin } 3800598e1a2SLisandro Dalcin ierr = PetscFree((*entities));CHKERRQ(ierr); 3810598e1a2SLisandro Dalcin PetscFunctionReturn(0); 3820598e1a2SLisandro Dalcin } 3830598e1a2SLisandro Dalcin 384730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 385730356f6SLisandro Dalcin { 386730356f6SLisandro Dalcin PetscErrorCode ierr; 3870598e1a2SLisandro Dalcin 388730356f6SLisandro Dalcin PetscFunctionBegin; 389730356f6SLisandro Dalcin ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr); 390730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 391730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 392730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 393730356f6SLisandro Dalcin PetscFunctionReturn(0); 394730356f6SLisandro Dalcin } 395730356f6SLisandro Dalcin 396730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 397730356f6SLisandro Dalcin { 398730356f6SLisandro Dalcin PetscInt index; 399730356f6SLisandro Dalcin PetscErrorCode ierr; 400730356f6SLisandro Dalcin 401730356f6SLisandro Dalcin PetscFunctionBegin; 402730356f6SLisandro Dalcin ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr); 403730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 404730356f6SLisandro Dalcin PetscFunctionReturn(0); 405730356f6SLisandro Dalcin } 406730356f6SLisandro Dalcin 4070598e1a2SLisandro Dalcin typedef struct { 4080598e1a2SLisandro Dalcin PetscInt *id; /* Node IDs */ 4090598e1a2SLisandro Dalcin double *xyz; /* Coordinates */ 4100598e1a2SLisandro Dalcin } GmshNodes; 4110598e1a2SLisandro Dalcin 4120598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes) 413730356f6SLisandro Dalcin { 414730356f6SLisandro Dalcin PetscErrorCode ierr; 415730356f6SLisandro Dalcin 416730356f6SLisandro Dalcin PetscFunctionBegin; 4170598e1a2SLisandro Dalcin ierr = PetscNew(nodes);CHKERRQ(ierr); 4180598e1a2SLisandro Dalcin ierr = PetscMalloc1(count*1, &(*nodes)->id);CHKERRQ(ierr); 4190598e1a2SLisandro Dalcin ierr = PetscMalloc1(count*3, &(*nodes)->xyz);CHKERRQ(ierr); 4200598e1a2SLisandro Dalcin PetscFunctionReturn(0); 421730356f6SLisandro Dalcin } 4220598e1a2SLisandro Dalcin 4230598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes) 4240598e1a2SLisandro Dalcin { 4250598e1a2SLisandro Dalcin PetscErrorCode ierr; 4260598e1a2SLisandro Dalcin PetscFunctionBegin; 4270598e1a2SLisandro Dalcin if (!*nodes) PetscFunctionReturn(0); 4280598e1a2SLisandro Dalcin ierr = PetscFree((*nodes)->id);CHKERRQ(ierr); 4290598e1a2SLisandro Dalcin ierr = PetscFree((*nodes)->xyz);CHKERRQ(ierr); 4300598e1a2SLisandro Dalcin ierr = PetscFree((*nodes));CHKERRQ(ierr); 431730356f6SLisandro Dalcin PetscFunctionReturn(0); 432730356f6SLisandro Dalcin } 433730356f6SLisandro Dalcin 434730356f6SLisandro Dalcin typedef struct { 4350598e1a2SLisandro Dalcin PetscInt id; /* Element ID */ 4360598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 437de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 4380598e1a2SLisandro Dalcin PetscInt numVerts; /* Size of vertex array */ 439de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 4400598e1a2SLisandro Dalcin PetscInt *nodes; /* Vertex/Node array */ 441de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 442de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 443de78e4feSLisandro Dalcin } GmshElement; 444de78e4feSLisandro Dalcin 4450598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements) 4460598e1a2SLisandro Dalcin { 4470598e1a2SLisandro Dalcin PetscErrorCode ierr; 4480598e1a2SLisandro Dalcin 4490598e1a2SLisandro Dalcin PetscFunctionBegin; 4500598e1a2SLisandro Dalcin ierr = PetscCalloc1(count, elements);CHKERRQ(ierr); 4510598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4520598e1a2SLisandro Dalcin } 4530598e1a2SLisandro Dalcin 4540598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsDestroy(GmshElement **elements) 4550598e1a2SLisandro Dalcin { 4560598e1a2SLisandro Dalcin PetscErrorCode ierr; 4570598e1a2SLisandro Dalcin 4580598e1a2SLisandro Dalcin PetscFunctionBegin; 4590598e1a2SLisandro Dalcin if (!*elements) PetscFunctionReturn(0); 4600598e1a2SLisandro Dalcin ierr = PetscFree(*elements);CHKERRQ(ierr); 4610598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4620598e1a2SLisandro Dalcin } 4630598e1a2SLisandro Dalcin 4640598e1a2SLisandro Dalcin typedef struct { 4650598e1a2SLisandro Dalcin PetscInt dim; 466066ea43fSLisandro Dalcin PetscInt order; 4670598e1a2SLisandro Dalcin GmshEntities *entities; 4680598e1a2SLisandro Dalcin PetscInt numNodes; 4690598e1a2SLisandro Dalcin GmshNodes *nodelist; 4700598e1a2SLisandro Dalcin PetscInt numElems; 4710598e1a2SLisandro Dalcin GmshElement *elements; 4720598e1a2SLisandro Dalcin PetscInt numVerts; 4730598e1a2SLisandro Dalcin PetscInt numCells; 4740598e1a2SLisandro Dalcin PetscInt *periodMap; 4750598e1a2SLisandro Dalcin PetscInt *vertexMap; 4760598e1a2SLisandro Dalcin PetscSegBuffer segbuf; 477*a45dabc8SMatthew G. Knepley PetscInt numRegions; 478*a45dabc8SMatthew G. Knepley PetscInt *regionTags; 479*a45dabc8SMatthew G. Knepley char **regionNames; 4800598e1a2SLisandro Dalcin } GmshMesh; 4810598e1a2SLisandro Dalcin 4820598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshCreate(GmshMesh **mesh) 4830598e1a2SLisandro Dalcin { 4840598e1a2SLisandro Dalcin PetscErrorCode ierr; 4850598e1a2SLisandro Dalcin 4860598e1a2SLisandro Dalcin PetscFunctionBegin; 4870598e1a2SLisandro Dalcin ierr = PetscNew(mesh);CHKERRQ(ierr); 4880598e1a2SLisandro Dalcin ierr = PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf);CHKERRQ(ierr); 4890598e1a2SLisandro Dalcin PetscFunctionReturn(0); 4900598e1a2SLisandro Dalcin } 4910598e1a2SLisandro Dalcin 4920598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh) 4930598e1a2SLisandro Dalcin { 494*a45dabc8SMatthew G. Knepley PetscInt r; 4950598e1a2SLisandro Dalcin PetscErrorCode ierr; 4960598e1a2SLisandro Dalcin 4970598e1a2SLisandro Dalcin PetscFunctionBegin; 4980598e1a2SLisandro Dalcin if (!*mesh) PetscFunctionReturn(0); 4990598e1a2SLisandro Dalcin ierr = GmshEntitiesDestroy(&(*mesh)->entities);CHKERRQ(ierr); 5000598e1a2SLisandro Dalcin ierr = GmshNodesDestroy(&(*mesh)->nodelist);CHKERRQ(ierr); 5010598e1a2SLisandro Dalcin ierr = GmshElementsDestroy(&(*mesh)->elements);CHKERRQ(ierr); 5020598e1a2SLisandro Dalcin ierr = PetscFree((*mesh)->periodMap);CHKERRQ(ierr); 5030598e1a2SLisandro Dalcin ierr = PetscFree((*mesh)->vertexMap);CHKERRQ(ierr); 5040598e1a2SLisandro Dalcin ierr = PetscSegBufferDestroy(&(*mesh)->segbuf);CHKERRQ(ierr); 505*a45dabc8SMatthew G. Knepley for (r = 0; r < (*mesh)->numRegions; ++r) {ierr = PetscFree((*mesh)->regionNames[r]);CHKERRQ(ierr);} 506*a45dabc8SMatthew G. Knepley ierr = PetscFree2((*mesh)->regionTags, (*mesh)->regionNames);CHKERRQ(ierr); 5070598e1a2SLisandro Dalcin ierr = PetscFree((*mesh));CHKERRQ(ierr); 5080598e1a2SLisandro Dalcin PetscFunctionReturn(0); 5090598e1a2SLisandro Dalcin } 5100598e1a2SLisandro Dalcin 5110598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh) 512de78e4feSLisandro Dalcin { 513698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 514698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 515de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 5160598e1a2SLisandro Dalcin int n, num, nid, snum; 5170598e1a2SLisandro Dalcin GmshNodes *nodes; 518de78e4feSLisandro Dalcin PetscErrorCode ierr; 519de78e4feSLisandro Dalcin 520de78e4feSLisandro Dalcin PetscFunctionBegin; 521de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 522698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 5230598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 5240598e1a2SLisandro Dalcin ierr = GmshNodesCreate(num, &nodes);CHKERRQ(ierr); 5250598e1a2SLisandro Dalcin mesh->numNodes = num; 5260598e1a2SLisandro Dalcin mesh->nodelist = nodes; 5270598e1a2SLisandro Dalcin for (n = 0; n < num; ++n) { 5280598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 52911c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 53011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 5310598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 53211c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 5330598e1a2SLisandro Dalcin nodes->id[n] = nid; 53411c56cb0SLisandro Dalcin } 53511c56cb0SLisandro Dalcin PetscFunctionReturn(0); 53611c56cb0SLisandro Dalcin } 53711c56cb0SLisandro Dalcin 538de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 539de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 540de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 541de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 5420598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v22(GmshFile* gmsh, GmshMesh *mesh) 54311c56cb0SLisandro Dalcin { 544698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 545698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 546698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 547de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 5480598e1a2SLisandro Dalcin int i, c, p, num, ibuf[1+4+1000], snum; 5490598e1a2SLisandro Dalcin int cellType, numElem, numVerts, numNodes, numTags; 550df032b82SMatthew G. Knepley GmshElement *elements; 5510598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 552df032b82SMatthew G. Knepley PetscErrorCode ierr; 553df032b82SMatthew G. Knepley 554df032b82SMatthew G. Knepley PetscFunctionBegin; 555feb237baSPierre Jolivet ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 556698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 5570598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 5580598e1a2SLisandro Dalcin ierr = GmshElementsCreate(num, &elements);CHKERRQ(ierr); 5590598e1a2SLisandro Dalcin mesh->numElems = num; 5600598e1a2SLisandro Dalcin mesh->elements = elements; 561698a79b9SLisandro Dalcin for (c = 0; c < num;) { 56211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 56311c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 5640598e1a2SLisandro Dalcin 5650598e1a2SLisandro Dalcin cellType = binary ? ibuf[0] : ibuf[1]; 5660598e1a2SLisandro Dalcin numElem = binary ? ibuf[1] : 1; 567df032b82SMatthew G. Knepley numTags = ibuf[2]; 5680598e1a2SLisandro Dalcin 5690598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 5700598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 5710598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 5720598e1a2SLisandro Dalcin 573df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 5740598e1a2SLisandro Dalcin GmshElement *element = elements + c; 5750598e1a2SLisandro Dalcin const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off; 5760598e1a2SLisandro Dalcin const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1; 5770598e1a2SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf+off, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 5780598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf+off, PETSC_ENUM, nint);CHKERRQ(ierr);} 5790598e1a2SLisandro Dalcin element->id = id[0]; 5800598e1a2SLisandro Dalcin element->dim = GmshCellMap[cellType].dim; 5810598e1a2SLisandro Dalcin element->cellType = cellType; 5820598e1a2SLisandro Dalcin element->numVerts = numVerts; 5830598e1a2SLisandro Dalcin element->numNodes = numNodes; 5840598e1a2SLisandro Dalcin element->numTags = PetscMin(numTags, 4); 5850598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 5860598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 5870598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = tags[p]; 588df032b82SMatthew G. Knepley } 589df032b82SMatthew G. Knepley } 590df032b82SMatthew G. Knepley PetscFunctionReturn(0); 591df032b82SMatthew G. Knepley } 5927d282ae0SMichael Lange 593de78e4feSLisandro Dalcin /* 594de78e4feSLisandro Dalcin $Entities 595de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 596de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 597de78e4feSLisandro Dalcin // points 598de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 599de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 600de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 601de78e4feSLisandro Dalcin ... 602de78e4feSLisandro Dalcin // curves 603de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 604de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 605de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 606de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 607de78e4feSLisandro Dalcin ... 608de78e4feSLisandro Dalcin // surfaces 609de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 610de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 611de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 612de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 613de78e4feSLisandro Dalcin ... 614de78e4feSLisandro Dalcin // volumes 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 numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 619de78e4feSLisandro Dalcin ... 620de78e4feSLisandro Dalcin $EndEntities 621de78e4feSLisandro Dalcin */ 6220598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh) 623de78e4feSLisandro Dalcin { 624698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 625698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 626698a79b9SLisandro Dalcin long index, num, lbuf[4]; 627730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 628698a79b9SLisandro Dalcin PetscInt count[4], i; 629a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 630de78e4feSLisandro Dalcin PetscErrorCode ierr; 631de78e4feSLisandro Dalcin 632de78e4feSLisandro Dalcin PetscFunctionBegin; 633698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 634698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 635698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 6360598e1a2SLisandro Dalcin ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr); 637de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 638730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 639730356f6SLisandro Dalcin ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 640730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 6410598e1a2SLisandro Dalcin ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 642de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 643de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 644de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 645de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 646698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 647de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 648730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 649de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 650de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 651de78e4feSLisandro Dalcin if (dim == 0) continue; 652de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 653de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 654698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 655de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 656730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 657de78e4feSLisandro Dalcin } 658de78e4feSLisandro Dalcin } 659de78e4feSLisandro Dalcin PetscFunctionReturn(0); 660de78e4feSLisandro Dalcin } 661de78e4feSLisandro Dalcin 662de78e4feSLisandro Dalcin /* 663de78e4feSLisandro Dalcin $Nodes 664de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 665de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 666de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 667de78e4feSLisandro Dalcin ... 668de78e4feSLisandro Dalcin ... 669de78e4feSLisandro Dalcin $EndNodes 670de78e4feSLisandro Dalcin */ 6710598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh) 672de78e4feSLisandro Dalcin { 673698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 674698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 6750598e1a2SLisandro Dalcin long block, node, n, numEntityBlocks, numTotalNodes, numNodes; 676de78e4feSLisandro Dalcin int info[3], nid; 6770598e1a2SLisandro Dalcin GmshNodes *nodes; 678de78e4feSLisandro Dalcin PetscErrorCode ierr; 679de78e4feSLisandro Dalcin 680de78e4feSLisandro Dalcin PetscFunctionBegin; 681de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 682de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 683de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 684de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 6850598e1a2SLisandro Dalcin ierr = GmshNodesCreate(numTotalNodes, &nodes);CHKERRQ(ierr); 6860598e1a2SLisandro Dalcin mesh->numNodes = numTotalNodes; 6870598e1a2SLisandro Dalcin mesh->nodelist = nodes; 6880598e1a2SLisandro Dalcin for (n = 0, block = 0; block < numEntityBlocks; ++block) { 689de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 690de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 691de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 692698a79b9SLisandro Dalcin if (gmsh->binary) { 693277f51e8SBarry Smith size_t nbytes = sizeof(int) + 3*sizeof(double); 694277f51e8SBarry Smith char *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */ 695698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 696de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 6970598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 698de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 6990598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 70030815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);} 70130815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 702de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 703de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 704de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 705de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 7060598e1a2SLisandro Dalcin nodes->id[n] = nid; 707de78e4feSLisandro Dalcin } 708de78e4feSLisandro Dalcin } else { 7090598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 7100598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 711de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 712de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 7130598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 714de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 7150598e1a2SLisandro Dalcin nodes->id[n] = nid; 716de78e4feSLisandro Dalcin } 717de78e4feSLisandro Dalcin } 718de78e4feSLisandro Dalcin } 719de78e4feSLisandro Dalcin PetscFunctionReturn(0); 720de78e4feSLisandro Dalcin } 721de78e4feSLisandro Dalcin 722de78e4feSLisandro Dalcin /* 723de78e4feSLisandro Dalcin $Elements 724de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 725de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 726de78e4feSLisandro Dalcin tag(int) numVert[...](int) 727de78e4feSLisandro Dalcin ... 728de78e4feSLisandro Dalcin ... 729de78e4feSLisandro Dalcin $EndElements 730de78e4feSLisandro Dalcin */ 7310598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh) 732de78e4feSLisandro Dalcin { 733698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 734698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 735de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 736de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 7370598e1a2SLisandro Dalcin int eid, dim, cellType, numVerts, numNodes, numTags; 738a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 739de78e4feSLisandro Dalcin GmshElement *elements; 7400598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 741de78e4feSLisandro Dalcin PetscErrorCode ierr; 742de78e4feSLisandro Dalcin 743de78e4feSLisandro Dalcin PetscFunctionBegin; 744de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 745de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 746de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 747de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 7480598e1a2SLisandro Dalcin ierr = GmshElementsCreate(numTotalElements, &elements);CHKERRQ(ierr); 7490598e1a2SLisandro Dalcin mesh->numElems = numTotalElements; 7500598e1a2SLisandro Dalcin mesh->elements = elements; 751de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 752de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 753de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 754de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 7550598e1a2SLisandro Dalcin ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr); 7560598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 7570598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 7580598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 759730356f6SLisandro Dalcin numTags = entity->numTags; 760de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 761de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 762698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 763de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 764de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 765de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 766de78e4feSLisandro Dalcin GmshElement *element = elements + c; 7670598e1a2SLisandro Dalcin const int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 7680598e1a2SLisandro Dalcin element->id = id[0]; 769de78e4feSLisandro Dalcin element->dim = dim; 770de78e4feSLisandro Dalcin element->cellType = cellType; 7710598e1a2SLisandro Dalcin element->numVerts = numVerts; 772de78e4feSLisandro Dalcin element->numNodes = numNodes; 773de78e4feSLisandro Dalcin element->numTags = numTags; 7740598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 7750598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 7760598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 777de78e4feSLisandro Dalcin } 778de78e4feSLisandro Dalcin } 779698a79b9SLisandro Dalcin PetscFunctionReturn(0); 780698a79b9SLisandro Dalcin } 781698a79b9SLisandro Dalcin 7820598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[]) 783698a79b9SLisandro Dalcin { 784698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 785698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 786698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 787698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 788698a79b9SLisandro Dalcin int numPeriodic, snum, i; 789698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 7900598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 791698a79b9SLisandro Dalcin PetscErrorCode ierr; 792698a79b9SLisandro Dalcin 793698a79b9SLisandro Dalcin PetscFunctionBegin; 794698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 795698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 796698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 7970598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 798698a79b9SLisandro Dalcin } else { 799698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 800698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 801698a79b9SLisandro Dalcin } 802698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 8039dddd249SSatish Balay int ibuf[3], correspondingDim = -1, correspondingTag = -1, primaryTag = -1, correspondingNode, primaryNode; 804698a79b9SLisandro Dalcin long j, nNodes; 805698a79b9SLisandro Dalcin double affine[16]; 806698a79b9SLisandro Dalcin 807698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 808698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 8099dddd249SSatish Balay snum = sscanf(line, "%d %d %d", &correspondingDim, &correspondingTag, &primaryTag); 8100598e1a2SLisandro Dalcin if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 811698a79b9SLisandro Dalcin } else { 812698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 813698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 8149dddd249SSatish Balay correspondingDim = ibuf[0]; correspondingTag = ibuf[1]; primaryTag = ibuf[2]; 815698a79b9SLisandro Dalcin } 8169dddd249SSatish Balay (void)correspondingDim; (void)correspondingTag; (void)primaryTag; /* unused */ 817698a79b9SLisandro Dalcin 818698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 819698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 820698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 8214c500f23SPierre Jolivet if (snum != 1) { /* discard transformation and try again */ 822698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 823698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 824698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 8250598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 826698a79b9SLisandro Dalcin } 827698a79b9SLisandro Dalcin } else { 828698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 829698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 8304c500f23SPierre Jolivet if (nNodes == -1) { /* discard transformation and try again */ 831698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 832698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 833698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 834698a79b9SLisandro Dalcin } 835698a79b9SLisandro Dalcin } 836698a79b9SLisandro Dalcin 837698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 838698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 839698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 8409dddd249SSatish Balay snum = sscanf(line, "%d %d", &correspondingNode, &primaryNode); 8410598e1a2SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 842698a79b9SLisandro Dalcin } else { 843698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 844698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 8459dddd249SSatish Balay correspondingNode = ibuf[0]; primaryNode = ibuf[1]; 846698a79b9SLisandro Dalcin } 8479dddd249SSatish Balay correspondingNode = (int) nodeMap[correspondingNode]; 8489dddd249SSatish Balay primaryNode = (int) nodeMap[primaryNode]; 8499dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 850698a79b9SLisandro Dalcin } 851698a79b9SLisandro Dalcin } 852698a79b9SLisandro Dalcin PetscFunctionReturn(0); 853698a79b9SLisandro Dalcin } 854698a79b9SLisandro Dalcin 8550598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 856698a79b9SLisandro Dalcin $Entities 857698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 858698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 859698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 860698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 861698a79b9SLisandro Dalcin ... 862698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 863698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 864698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 865698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 866698a79b9SLisandro Dalcin ... 867698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 868698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 869698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 870698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 871698a79b9SLisandro Dalcin ... 872698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 873698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 874698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 875698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 876698a79b9SLisandro Dalcin ... 877698a79b9SLisandro Dalcin $EndEntities 878698a79b9SLisandro Dalcin */ 8790598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh) 880698a79b9SLisandro Dalcin { 881698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 882698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 883698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 884698a79b9SLisandro Dalcin PetscErrorCode ierr; 885698a79b9SLisandro Dalcin 886698a79b9SLisandro Dalcin PetscFunctionBegin; 887698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr); 8880598e1a2SLisandro Dalcin ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr); 889698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 890698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 891698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr); 8920598e1a2SLisandro Dalcin ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 893698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr); 894698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 895698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 896698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 897698a79b9SLisandro Dalcin entity->numTags = PetscMin(numTags, 4); 898698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 899698a79b9SLisandro Dalcin if (dim == 0) continue; 900698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 901698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 902698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 903698a79b9SLisandro Dalcin } 904698a79b9SLisandro Dalcin } 905698a79b9SLisandro Dalcin PetscFunctionReturn(0); 906698a79b9SLisandro Dalcin } 907698a79b9SLisandro Dalcin 90833a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 909698a79b9SLisandro Dalcin $Nodes 910698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 911698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 912698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 913698a79b9SLisandro Dalcin nodeTag(size_t) 914698a79b9SLisandro Dalcin ... 915698a79b9SLisandro Dalcin x(double) y(double) z(double) 916698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 917698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 918698a79b9SLisandro Dalcin ... 919698a79b9SLisandro Dalcin ... 920698a79b9SLisandro Dalcin $EndNodes 921698a79b9SLisandro Dalcin */ 9220598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh) 923698a79b9SLisandro Dalcin { 924698a79b9SLisandro Dalcin int info[3]; 9250598e1a2SLisandro Dalcin PetscInt sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, block, node; 9260598e1a2SLisandro Dalcin GmshNodes *nodes; 927698a79b9SLisandro Dalcin PetscErrorCode ierr; 928698a79b9SLisandro Dalcin 929698a79b9SLisandro Dalcin PetscFunctionBegin; 930698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 9310598e1a2SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 9320598e1a2SLisandro Dalcin ierr = GmshNodesCreate(numNodes, &nodes);CHKERRQ(ierr); 9330598e1a2SLisandro Dalcin mesh->numNodes = numNodes; 9340598e1a2SLisandro Dalcin mesh->nodelist = nodes; 935698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 936698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 937698a79b9SLisandro Dalcin if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 9380598e1a2SLisandro Dalcin ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr); 9390598e1a2SLisandro Dalcin ierr = GmshReadSize(gmsh, nodes->id+node, numNodesBlock);CHKERRQ(ierr); 9400598e1a2SLisandro Dalcin ierr = GmshReadDouble(gmsh, nodes->xyz+node*3, numNodesBlock*3);CHKERRQ(ierr); 941698a79b9SLisandro Dalcin } 9420598e1a2SLisandro Dalcin gmsh->nodeStart = sizes[2]; 9430598e1a2SLisandro Dalcin gmsh->nodeEnd = sizes[3]+1; 944698a79b9SLisandro Dalcin PetscFunctionReturn(0); 945698a79b9SLisandro Dalcin } 946698a79b9SLisandro Dalcin 94733a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 948698a79b9SLisandro Dalcin $Elements 949698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 950698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 951698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 952698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 953698a79b9SLisandro Dalcin ... 954698a79b9SLisandro Dalcin ... 955698a79b9SLisandro Dalcin $EndElements 956698a79b9SLisandro Dalcin */ 9570598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh) 958698a79b9SLisandro Dalcin { 9590598e1a2SLisandro Dalcin int info[3], eid, dim, cellType; 9600598e1a2SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p; 961698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 962698a79b9SLisandro Dalcin GmshElement *elements; 9630598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 964698a79b9SLisandro Dalcin PetscErrorCode ierr; 965698a79b9SLisandro Dalcin 966698a79b9SLisandro Dalcin PetscFunctionBegin; 967698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 968698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 9690598e1a2SLisandro Dalcin ierr = GmshElementsCreate(numElements, &elements);CHKERRQ(ierr); 9700598e1a2SLisandro Dalcin mesh->numElems = numElements; 9710598e1a2SLisandro Dalcin mesh->elements = elements; 972698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 973698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 974698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 9750598e1a2SLisandro Dalcin ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr); 9760598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 9770598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 9780598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 979698a79b9SLisandro Dalcin numTags = entity->numTags; 980698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr); 981698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr); 982698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr); 983698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 984698a79b9SLisandro Dalcin GmshElement *element = elements + c; 9850598e1a2SLisandro Dalcin const PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 986698a79b9SLisandro Dalcin element->id = id[0]; 987698a79b9SLisandro Dalcin element->dim = dim; 988698a79b9SLisandro Dalcin element->cellType = cellType; 9890598e1a2SLisandro Dalcin element->numVerts = numVerts; 990698a79b9SLisandro Dalcin element->numNodes = numNodes; 991698a79b9SLisandro Dalcin element->numTags = numTags; 9920598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 9930598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 9940598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 995698a79b9SLisandro Dalcin } 996698a79b9SLisandro Dalcin } 997698a79b9SLisandro Dalcin PetscFunctionReturn(0); 998698a79b9SLisandro Dalcin } 999698a79b9SLisandro Dalcin 10000598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 1001698a79b9SLisandro Dalcin $Periodic 1002698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 10039dddd249SSatish Balay entityDim(int) entityTag(int) entityTagPrimary(int) 1004698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 1005698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 10069dddd249SSatish Balay nodeTag(size_t) nodeTagPrimary(size_t) 1007698a79b9SLisandro Dalcin ... 1008698a79b9SLisandro Dalcin ... 1009698a79b9SLisandro Dalcin $EndPeriodic 1010698a79b9SLisandro Dalcin */ 10110598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[]) 1012698a79b9SLisandro Dalcin { 1013698a79b9SLisandro Dalcin int info[3]; 1014698a79b9SLisandro Dalcin double dbuf[16]; 10150598e1a2SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 10160598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 1017698a79b9SLisandro Dalcin PetscErrorCode ierr; 1018698a79b9SLisandro Dalcin 1019698a79b9SLisandro Dalcin PetscFunctionBegin; 1020698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr); 1021698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 1022698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 1023698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr); 1024698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr); 1025698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr); 1026698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr); 1027698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr); 1028698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 10299dddd249SSatish Balay PetscInt correspondingNode = nodeMap[nodeTags[node*2+0]]; 10309dddd249SSatish Balay PetscInt primaryNode = nodeMap[nodeTags[node*2+1]]; 10319dddd249SSatish Balay periodicMap[correspondingNode] = primaryNode; 1032698a79b9SLisandro Dalcin } 1033698a79b9SLisandro Dalcin } 1034698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1035698a79b9SLisandro Dalcin } 1036698a79b9SLisandro Dalcin 10370598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 1038d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 1039d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 1040d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 1041d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 1042d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 1043d6689ca9SLisandro Dalcin $EndMeshFormat 1044d6689ca9SLisandro Dalcin */ 10450598e1a2SLisandro Dalcin static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh) 1046698a79b9SLisandro Dalcin { 1047698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1048698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 1049698a79b9SLisandro Dalcin float version; 1050698a79b9SLisandro Dalcin PetscErrorCode ierr; 1051698a79b9SLisandro Dalcin 1052698a79b9SLisandro Dalcin PetscFunctionBegin; 1053698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 1054698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 10550598e1a2SLisandro Dalcin if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 10560598e1a2SLisandro Dalcin if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 10570598e1a2SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 10580598e1a2SLisandro Dalcin if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 10590598e1a2SLisandro Dalcin if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII"); 10600598e1a2SLisandro Dalcin if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary"); 1061af7a0b9fSSatish Balay fileFormat = (int)roundf(version*10); 10620598e1a2SLisandro Dalcin if (fileFormat <= 40 && dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize); 10630598e1a2SLisandro Dalcin if (fileFormat >= 41 && dataSize != sizeof(int) && dataSize != sizeof(PetscInt64)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Data size %d is not valid for a Gmsh file", dataSize); 1064698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 1065698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 1066698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 1067698a79b9SLisandro Dalcin if (gmsh->binary) { 1068698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr); 1069698a79b9SLisandro Dalcin if (checkEndian != 1) { 1070698a79b9SLisandro Dalcin ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr); 10710598e1a2SLisandro Dalcin if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line); 1072698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 1073698a79b9SLisandro Dalcin } 1074698a79b9SLisandro Dalcin } 1075698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1076698a79b9SLisandro Dalcin } 1077698a79b9SLisandro Dalcin 10781f643af3SLisandro Dalcin /* 10791f643af3SLisandro Dalcin PhysicalNames 10801f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 10811f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 10821f643af3SLisandro Dalcin ... 10831f643af3SLisandro Dalcin $EndPhysicalNames 10841f643af3SLisandro Dalcin */ 1085*a45dabc8SMatthew G. Knepley static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh, GmshMesh *mesh) 1086698a79b9SLisandro Dalcin { 10871f643af3SLisandro Dalcin char line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q; 1088*a45dabc8SMatthew G. Knepley int snum, region, dim, tag; 1089698a79b9SLisandro Dalcin PetscErrorCode ierr; 1090698a79b9SLisandro Dalcin 1091698a79b9SLisandro Dalcin PetscFunctionBegin; 1092698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1093*a45dabc8SMatthew G. Knepley snum = sscanf(line, "%d", ®ion); 1094*a45dabc8SMatthew G. Knepley mesh->numRegions = region; 10950598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 1096*a45dabc8SMatthew G. Knepley ierr = PetscMalloc2(mesh->numRegions, &mesh->regionTags, mesh->numRegions, &mesh->regionNames);CHKERRQ(ierr); 1097*a45dabc8SMatthew G. Knepley for (region = 0; region < mesh->numRegions; ++region) { 10981f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 10991f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 11000598e1a2SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11011f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr); 11021f643af3SLisandro Dalcin ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr); 11030598e1a2SLisandro Dalcin if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11041f643af3SLisandro Dalcin ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr); 11050598e1a2SLisandro Dalcin if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 11061f643af3SLisandro Dalcin ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr); 1107*a45dabc8SMatthew G. Knepley mesh->regionTags[region] = tag; 1108*a45dabc8SMatthew G. Knepley ierr = PetscStrallocpy(name, &mesh->regionNames[region]);CHKERRQ(ierr); 1109698a79b9SLisandro Dalcin } 1110de78e4feSLisandro Dalcin PetscFunctionReturn(0); 1111de78e4feSLisandro Dalcin } 1112de78e4feSLisandro Dalcin 11130598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh) 111496ca5757SLisandro Dalcin { 11150598e1a2SLisandro Dalcin PetscErrorCode ierr; 11160598e1a2SLisandro Dalcin 11170598e1a2SLisandro Dalcin PetscFunctionBegin; 11180598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11190598e1a2SLisandro Dalcin case 41: ierr = GmshReadEntities_v41(gmsh, mesh);CHKERRQ(ierr); break; 11200598e1a2SLisandro Dalcin default: ierr = GmshReadEntities_v40(gmsh, mesh);CHKERRQ(ierr); break; 112196ca5757SLisandro Dalcin } 11220598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11230598e1a2SLisandro Dalcin } 11240598e1a2SLisandro Dalcin 11250598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh) 11260598e1a2SLisandro Dalcin { 11270598e1a2SLisandro Dalcin PetscErrorCode ierr; 11280598e1a2SLisandro Dalcin 11290598e1a2SLisandro Dalcin PetscFunctionBegin; 11300598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11310598e1a2SLisandro Dalcin case 41: ierr = GmshReadNodes_v41(gmsh, mesh);CHKERRQ(ierr); break; 11320598e1a2SLisandro Dalcin case 40: ierr = GmshReadNodes_v40(gmsh, mesh);CHKERRQ(ierr); break; 11330598e1a2SLisandro Dalcin default: ierr = GmshReadNodes_v22(gmsh, mesh);CHKERRQ(ierr); break; 11340598e1a2SLisandro Dalcin } 11350598e1a2SLisandro Dalcin 11360598e1a2SLisandro Dalcin { /* Gmsh v2.2/v4.0 does not provide min/max node tags */ 11370598e1a2SLisandro Dalcin if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) { 11380598e1a2SLisandro Dalcin PetscInt tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n; 11390598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11400598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11410598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 11420598e1a2SLisandro Dalcin tagMin = PetscMin(tag, tagMin); 11430598e1a2SLisandro Dalcin tagMax = PetscMax(tag, tagMax); 11440598e1a2SLisandro Dalcin } 11450598e1a2SLisandro Dalcin gmsh->nodeStart = tagMin; 11460598e1a2SLisandro Dalcin gmsh->nodeEnd = tagMax+1; 11470598e1a2SLisandro Dalcin } 11480598e1a2SLisandro Dalcin } 11490598e1a2SLisandro Dalcin 11500598e1a2SLisandro Dalcin { /* Support for sparse node tags */ 11510598e1a2SLisandro Dalcin PetscInt n, t; 11520598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 11530598e1a2SLisandro Dalcin ierr = PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf);CHKERRQ(ierr); 11540598e1a2SLisandro Dalcin for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT; 11550598e1a2SLisandro Dalcin gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart; 11560598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 11570598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 11580598e1a2SLisandro Dalcin if (gmsh->nodeMap[tag] >= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %D", tag); 11590598e1a2SLisandro Dalcin gmsh->nodeMap[tag] = n; 11600598e1a2SLisandro Dalcin } 11610598e1a2SLisandro Dalcin } 11620598e1a2SLisandro Dalcin PetscFunctionReturn(0); 11630598e1a2SLisandro Dalcin } 11640598e1a2SLisandro Dalcin 11650598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh) 11660598e1a2SLisandro Dalcin { 11670598e1a2SLisandro Dalcin PetscErrorCode ierr; 11680598e1a2SLisandro Dalcin 11690598e1a2SLisandro Dalcin PetscFunctionBegin; 11700598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 11710598e1a2SLisandro Dalcin case 41: ierr = GmshReadElements_v41(gmsh, mesh);CHKERRQ(ierr); break; 11720598e1a2SLisandro Dalcin case 40: ierr = GmshReadElements_v40(gmsh, mesh);CHKERRQ(ierr); break; 11730598e1a2SLisandro Dalcin default: ierr = GmshReadElements_v22(gmsh, mesh);CHKERRQ(ierr); break; 11740598e1a2SLisandro Dalcin } 11750598e1a2SLisandro Dalcin 11760598e1a2SLisandro Dalcin { /* Reorder elements by codimension and polytope type */ 11770598e1a2SLisandro Dalcin PetscInt ne = mesh->numElems; 11780598e1a2SLisandro Dalcin GmshElement *elements = mesh->elements; 1179066ea43fSLisandro Dalcin PetscInt keymap[GMSH_NUM_POLYTOPES], nk = 0; 1180066ea43fSLisandro Dalcin PetscInt offset[GMSH_NUM_POLYTOPES+1], e, k; 11810598e1a2SLisandro Dalcin 1182066ea43fSLisandro Dalcin for (k = 0; k < GMSH_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT; 11830598e1a2SLisandro Dalcin ierr = PetscMemzero(offset,sizeof(offset));CHKERRQ(ierr); 11840598e1a2SLisandro Dalcin 1185066ea43fSLisandro Dalcin keymap[GMSH_TET] = nk++; 1186066ea43fSLisandro Dalcin keymap[GMSH_HEX] = nk++; 1187066ea43fSLisandro Dalcin keymap[GMSH_PRI] = nk++; 1188066ea43fSLisandro Dalcin keymap[GMSH_PYR] = nk++; 1189066ea43fSLisandro Dalcin keymap[GMSH_TRI] = nk++; 1190066ea43fSLisandro Dalcin keymap[GMSH_QUA] = nk++; 1191066ea43fSLisandro Dalcin keymap[GMSH_SEG] = nk++; 1192066ea43fSLisandro Dalcin keymap[GMSH_VTX] = nk++; 11930598e1a2SLisandro Dalcin 11940598e1a2SLisandro Dalcin ierr = GmshElementsCreate(mesh->numElems, &mesh->elements);CHKERRQ(ierr); 11950598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope] 11960598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) offset[1+key(e)]++; 11970598e1a2SLisandro Dalcin for (k = 1; k < nk; ++k) offset[k] += offset[k-1]; 11980598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e]; 11990598e1a2SLisandro Dalcin #undef key 12000598e1a2SLisandro Dalcin ierr = GmshElementsDestroy(&elements);CHKERRQ(ierr); 1201066ea43fSLisandro Dalcin } 12020598e1a2SLisandro Dalcin 1203066ea43fSLisandro Dalcin { /* Mesh dimension and order */ 1204066ea43fSLisandro Dalcin GmshElement *elem = mesh->numElems ? mesh->elements : NULL; 1205066ea43fSLisandro Dalcin mesh->dim = elem ? GmshCellMap[elem->cellType].dim : 0; 1206066ea43fSLisandro Dalcin mesh->order = elem ? GmshCellMap[elem->cellType].order : 0; 12070598e1a2SLisandro Dalcin } 12080598e1a2SLisandro Dalcin 12090598e1a2SLisandro Dalcin { 12100598e1a2SLisandro Dalcin PetscBT vtx; 12110598e1a2SLisandro Dalcin PetscInt dim = mesh->dim, e, n, v; 12120598e1a2SLisandro Dalcin 12130598e1a2SLisandro Dalcin ierr = PetscBTCreate(mesh->numNodes, &vtx);CHKERRQ(ierr); 12140598e1a2SLisandro Dalcin 12150598e1a2SLisandro Dalcin /* Compute number of cells and set of vertices */ 12160598e1a2SLisandro Dalcin mesh->numCells = 0; 12170598e1a2SLisandro Dalcin for (e = 0; e < mesh->numElems; ++e) { 12180598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 12190598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) mesh->numCells++; 12200598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; v++) { 12210598e1a2SLisandro Dalcin ierr = PetscBTSet(vtx, elem->nodes[v]);CHKERRQ(ierr); 12220598e1a2SLisandro Dalcin } 12230598e1a2SLisandro Dalcin } 12240598e1a2SLisandro Dalcin 12250598e1a2SLisandro Dalcin /* Compute numbering for vertices */ 12260598e1a2SLisandro Dalcin mesh->numVerts = 0; 12270598e1a2SLisandro Dalcin ierr = PetscMalloc1(mesh->numNodes, &mesh->vertexMap);CHKERRQ(ierr); 12280598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12290598e1a2SLisandro Dalcin mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT; 12300598e1a2SLisandro Dalcin 12310598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&vtx);CHKERRQ(ierr); 12320598e1a2SLisandro Dalcin } 12330598e1a2SLisandro Dalcin PetscFunctionReturn(0); 12340598e1a2SLisandro Dalcin } 12350598e1a2SLisandro Dalcin 12360598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh) 12370598e1a2SLisandro Dalcin { 12380598e1a2SLisandro Dalcin PetscInt n; 12390598e1a2SLisandro Dalcin PetscErrorCode ierr; 12400598e1a2SLisandro Dalcin 12410598e1a2SLisandro Dalcin PetscFunctionBegin; 12420598e1a2SLisandro Dalcin ierr = PetscMalloc1(mesh->numNodes, &mesh->periodMap);CHKERRQ(ierr); 12430598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n; 12440598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 12450598e1a2SLisandro Dalcin case 41: ierr = GmshReadPeriodic_v41(gmsh, mesh->periodMap);CHKERRQ(ierr); break; 12460598e1a2SLisandro Dalcin default: ierr = GmshReadPeriodic_v40(gmsh, mesh->periodMap);CHKERRQ(ierr); break; 12470598e1a2SLisandro Dalcin } 12480598e1a2SLisandro Dalcin 12499dddd249SSatish Balay /* Find canonical primary nodes */ 12500598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12510598e1a2SLisandro Dalcin while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]]) 12520598e1a2SLisandro Dalcin mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]]; 12530598e1a2SLisandro Dalcin 12549dddd249SSatish Balay /* Renumber vertices (filter out correspondings) */ 12550598e1a2SLisandro Dalcin mesh->numVerts = 0; 12560598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12570598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12589dddd249SSatish Balay if (mesh->periodMap[n] == n) /* is primary */ 12590598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->numVerts++; 12600598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 12610598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 12629dddd249SSatish Balay if (mesh->periodMap[n] != n) /* is corresponding */ 12630598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]]; 12640598e1a2SLisandro Dalcin PetscFunctionReturn(0); 12650598e1a2SLisandro Dalcin } 12660598e1a2SLisandro Dalcin 1267066ea43fSLisandro Dalcin #define DM_POLYTOPE_VERTEX DM_POLYTOPE_POINT 1268066ea43fSLisandro Dalcin #define DM_POLYTOPE_PYRAMID DM_POLYTOPE_UNKNOWN 1269066ea43fSLisandro Dalcin static const DMPolytopeType DMPolytopeMap[] = { 1270066ea43fSLisandro Dalcin /* GMSH_VTX */ DM_POLYTOPE_VERTEX, 1271066ea43fSLisandro Dalcin /* GMSH_SEG */ DM_POLYTOPE_SEGMENT, 1272066ea43fSLisandro Dalcin /* GMSH_TRI */ DM_POLYTOPE_TRIANGLE, 1273066ea43fSLisandro Dalcin /* GMSH_QUA */ DM_POLYTOPE_QUADRILATERAL, 1274066ea43fSLisandro Dalcin /* GMSH_TET */ DM_POLYTOPE_TETRAHEDRON, 1275066ea43fSLisandro Dalcin /* GMSH_HEX */ DM_POLYTOPE_HEXAHEDRON, 1276066ea43fSLisandro Dalcin /* GMSH_PRI */ DM_POLYTOPE_TRI_PRISM, 1277066ea43fSLisandro Dalcin /* GMSH_PYR */ DM_POLYTOPE_PYRAMID, 1278066ea43fSLisandro Dalcin DM_POLYTOPE_UNKNOWN 1279066ea43fSLisandro Dalcin }; 12800598e1a2SLisandro Dalcin 12810598e1a2SLisandro Dalcin PETSC_STATIC_INLINE DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType) 12820598e1a2SLisandro Dalcin { 1283066ea43fSLisandro Dalcin return DMPolytopeMap[GmshCellMap[cellType].polytope]; 1284066ea43fSLisandro Dalcin } 1285066ea43fSLisandro Dalcin 1286066ea43fSLisandro Dalcin static PetscErrorCode GmshCreateFE(MPI_Comm comm, const char prefix[], PetscBool isSimplex, PetscBool continuity, PetscDTNodeType nodeType, PetscInt dim, PetscInt Nc, PetscInt k, PetscFE *fem) 1287066ea43fSLisandro Dalcin { 1288066ea43fSLisandro Dalcin DM K; 1289066ea43fSLisandro Dalcin PetscSpace P; 1290066ea43fSLisandro Dalcin PetscDualSpace Q; 1291066ea43fSLisandro Dalcin PetscQuadrature q, fq; 1292066ea43fSLisandro Dalcin PetscBool isTensor = isSimplex ? PETSC_FALSE : PETSC_TRUE; 1293066ea43fSLisandro Dalcin PetscBool endpoint = PETSC_TRUE; 1294066ea43fSLisandro Dalcin char name[32]; 1295066ea43fSLisandro Dalcin PetscErrorCode ierr; 1296066ea43fSLisandro Dalcin 1297066ea43fSLisandro Dalcin PetscFunctionBegin; 1298066ea43fSLisandro Dalcin /* Create space */ 1299066ea43fSLisandro Dalcin ierr = PetscSpaceCreate(comm, &P);CHKERRQ(ierr); 1300066ea43fSLisandro Dalcin ierr = PetscSpaceSetType(P, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 1301066ea43fSLisandro Dalcin ierr = PetscSpacePolynomialSetTensor(P, isTensor);CHKERRQ(ierr); 1302066ea43fSLisandro Dalcin ierr = PetscSpaceSetNumComponents(P, Nc);CHKERRQ(ierr); 1303066ea43fSLisandro Dalcin ierr = PetscSpaceSetNumVariables(P, dim);CHKERRQ(ierr); 1304066ea43fSLisandro Dalcin ierr = PetscSpaceSetDegree(P, k, PETSC_DETERMINE);CHKERRQ(ierr); 1305066ea43fSLisandro Dalcin if (prefix) { 1306066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) P, prefix);CHKERRQ(ierr); 1307066ea43fSLisandro Dalcin ierr = PetscSpaceSetFromOptions(P);CHKERRQ(ierr); 1308066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) P, NULL);CHKERRQ(ierr); 1309066ea43fSLisandro Dalcin ierr = PetscSpaceGetDegree(P, &k, NULL);CHKERRQ(ierr); 1310066ea43fSLisandro Dalcin } 1311066ea43fSLisandro Dalcin ierr = PetscSpaceSetUp(P);CHKERRQ(ierr); 1312066ea43fSLisandro Dalcin /* Create dual space */ 1313066ea43fSLisandro Dalcin ierr = PetscDualSpaceCreate(comm, &Q);CHKERRQ(ierr); 1314066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetType(Q, PETSCDUALSPACELAGRANGE);CHKERRQ(ierr); 1315066ea43fSLisandro Dalcin ierr = PetscDualSpaceLagrangeSetTensor(Q, isTensor);CHKERRQ(ierr); 1316066ea43fSLisandro Dalcin ierr = PetscDualSpaceLagrangeSetContinuity(Q, continuity);CHKERRQ(ierr); 1317066ea43fSLisandro Dalcin ierr = PetscDualSpaceLagrangeSetNodeType(Q, nodeType, endpoint, 0);CHKERRQ(ierr); 1318066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetNumComponents(Q, Nc);CHKERRQ(ierr); 1319066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetOrder(Q, k);CHKERRQ(ierr); 1320066ea43fSLisandro Dalcin ierr = PetscDualSpaceCreateReferenceCell(Q, dim, isSimplex, &K);CHKERRQ(ierr); 1321066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetDM(Q, K);CHKERRQ(ierr); 1322066ea43fSLisandro Dalcin ierr = DMDestroy(&K);CHKERRQ(ierr); 1323066ea43fSLisandro Dalcin if (prefix) { 1324066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, prefix);CHKERRQ(ierr); 1325066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetFromOptions(Q);CHKERRQ(ierr); 1326066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) Q, NULL);CHKERRQ(ierr); 1327066ea43fSLisandro Dalcin } 1328066ea43fSLisandro Dalcin ierr = PetscDualSpaceSetUp(Q);CHKERRQ(ierr); 1329066ea43fSLisandro Dalcin /* Create quadrature */ 1330066ea43fSLisandro Dalcin if (isSimplex) { 1331066ea43fSLisandro Dalcin ierr = PetscDTStroudConicalQuadrature(dim, 1, k+1, -1, +1, &q);CHKERRQ(ierr); 1332066ea43fSLisandro Dalcin ierr = PetscDTStroudConicalQuadrature(dim-1, 1, k+1, -1, +1, &fq);CHKERRQ(ierr); 1333066ea43fSLisandro Dalcin } else { 1334066ea43fSLisandro Dalcin ierr = PetscDTGaussTensorQuadrature(dim, 1, k+1, -1, +1, &q);CHKERRQ(ierr); 1335066ea43fSLisandro Dalcin ierr = PetscDTGaussTensorQuadrature(dim-1, 1, k+1, -1, +1, &fq);CHKERRQ(ierr); 1336066ea43fSLisandro Dalcin } 1337066ea43fSLisandro Dalcin /* Create finite element */ 1338066ea43fSLisandro Dalcin ierr = PetscFECreate(comm, fem);CHKERRQ(ierr); 1339066ea43fSLisandro Dalcin ierr = PetscFESetType(*fem, PETSCFEBASIC);CHKERRQ(ierr); 1340066ea43fSLisandro Dalcin ierr = PetscFESetNumComponents(*fem, Nc);CHKERRQ(ierr); 1341066ea43fSLisandro Dalcin ierr = PetscFESetBasisSpace(*fem, P);CHKERRQ(ierr); 1342066ea43fSLisandro Dalcin ierr = PetscFESetDualSpace(*fem, Q);CHKERRQ(ierr); 1343066ea43fSLisandro Dalcin ierr = PetscFESetQuadrature(*fem, q);CHKERRQ(ierr); 1344066ea43fSLisandro Dalcin ierr = PetscFESetFaceQuadrature(*fem, fq);CHKERRQ(ierr); 1345066ea43fSLisandro Dalcin if (prefix) { 1346066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, prefix);CHKERRQ(ierr); 1347066ea43fSLisandro Dalcin ierr = PetscFESetFromOptions(*fem);CHKERRQ(ierr); 1348066ea43fSLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject) *fem, NULL);CHKERRQ(ierr); 1349066ea43fSLisandro Dalcin } 1350066ea43fSLisandro Dalcin ierr = PetscFESetUp(*fem);CHKERRQ(ierr); 1351066ea43fSLisandro Dalcin /* Cleanup */ 1352066ea43fSLisandro Dalcin ierr = PetscSpaceDestroy(&P);CHKERRQ(ierr); 1353066ea43fSLisandro Dalcin ierr = PetscDualSpaceDestroy(&Q);CHKERRQ(ierr); 1354066ea43fSLisandro Dalcin ierr = PetscQuadratureDestroy(&q);CHKERRQ(ierr); 1355066ea43fSLisandro Dalcin ierr = PetscQuadratureDestroy(&fq);CHKERRQ(ierr); 1356066ea43fSLisandro Dalcin /* Set finite element name */ 1357066ea43fSLisandro Dalcin ierr = PetscSNPrintf(name, sizeof(name), "%s%D", isSimplex? "P" : "Q", k);CHKERRQ(ierr); 1358066ea43fSLisandro Dalcin ierr = PetscFESetName(*fem, name);CHKERRQ(ierr); 1359066ea43fSLisandro Dalcin PetscFunctionReturn(0); 136096ca5757SLisandro Dalcin } 136196ca5757SLisandro Dalcin 1362d6689ca9SLisandro Dalcin /*@C 1363d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 1364d6689ca9SLisandro Dalcin 1365d6689ca9SLisandro Dalcin + comm - The MPI communicator 1366d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 1367d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 1368d6689ca9SLisandro Dalcin 1369d6689ca9SLisandro Dalcin Output Parameter: 1370d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 1371d6689ca9SLisandro Dalcin 1372d6689ca9SLisandro Dalcin Level: beginner 1373d6689ca9SLisandro Dalcin 1374d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 1375d6689ca9SLisandro Dalcin @*/ 1376d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1377d6689ca9SLisandro Dalcin { 1378d6689ca9SLisandro Dalcin PetscViewer viewer; 1379d6689ca9SLisandro Dalcin PetscMPIInt rank; 1380d6689ca9SLisandro Dalcin int fileType; 1381d6689ca9SLisandro Dalcin PetscViewerType vtype; 1382d6689ca9SLisandro Dalcin PetscErrorCode ierr; 1383d6689ca9SLisandro Dalcin 1384d6689ca9SLisandro Dalcin PetscFunctionBegin; 1385ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 1386d6689ca9SLisandro Dalcin 1387d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 1388dd400576SPatrick Sanan if (rank == 0) { 13890598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1390d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1391d6689ca9SLisandro Dalcin int snum; 1392d6689ca9SLisandro Dalcin float version; 1393d6689ca9SLisandro Dalcin 1394580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 1395d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr); 1396d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 1397d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr); 1398d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr); 1399d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 1400d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1401d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 1402d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 1403d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 14040598e1a2SLisandro Dalcin if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 14050598e1a2SLisandro Dalcin if (version < 2.2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at least 2.2", (double)version); 14060598e1a2SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 14070598e1a2SLisandro Dalcin if (version > 4.1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f must be at most 4.1", (double)version); 1408d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr); 1409d6689ca9SLisandro Dalcin } 1410ffc4695bSBarry Smith ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 1411d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 1412d6689ca9SLisandro Dalcin 1413d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 1414d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 1415d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 1416d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 1417d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 1418d6689ca9SLisandro Dalcin ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 1419d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1420d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 1421d6689ca9SLisandro Dalcin } 1422d6689ca9SLisandro Dalcin 1423331830f3SMatthew G. Knepley /*@ 14247d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 1425331830f3SMatthew G. Knepley 1426d083f849SBarry Smith Collective 1427331830f3SMatthew G. Knepley 1428331830f3SMatthew G. Knepley Input Parameters: 1429331830f3SMatthew G. Knepley + comm - The MPI communicator 1430331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1431331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1432331830f3SMatthew G. Knepley 1433331830f3SMatthew G. Knepley Output Parameter: 1434331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1435331830f3SMatthew G. Knepley 1436698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1437331830f3SMatthew G. Knepley 1438331830f3SMatthew G. Knepley Level: beginner 1439331830f3SMatthew G. Knepley 1440331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 1441331830f3SMatthew G. Knepley @*/ 1442331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1443331830f3SMatthew G. Knepley { 14440598e1a2SLisandro Dalcin GmshMesh *mesh = NULL; 144511c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 14460598e1a2SLisandro Dalcin PetscBT periodicVerts = NULL; 14470598e1a2SLisandro Dalcin PetscBT periodicCells = NULL; 1448066ea43fSLisandro Dalcin DM cdm; 1449331830f3SMatthew G. Knepley PetscSection coordSection; 1450331830f3SMatthew G. Knepley Vec coordinates; 1451*a45dabc8SMatthew G. Knepley DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL, marker = NULL, *regionSets; 1452066ea43fSLisandro Dalcin PetscInt dim = 0, coordDim = -1, order = 0; 14530598e1a2SLisandro Dalcin PetscInt numNodes = 0, numElems = 0, numVerts = 0, numCells = 0; 1454066ea43fSLisandro Dalcin PetscInt cell, cone[8], e, n, v, d; 1455*a45dabc8SMatthew G. Knepley PetscBool binary, usemarker = PETSC_FALSE, useregions = PETSC_FALSE; 14560598e1a2SLisandro Dalcin PetscBool hybrid = interpolate, periodic = PETSC_TRUE; 1457066ea43fSLisandro Dalcin PetscBool highOrder = PETSC_TRUE, highOrderSet, project = PETSC_FALSE; 1458066ea43fSLisandro Dalcin PetscBool isSimplex = PETSC_FALSE, isHybrid = PETSC_FALSE, hasTetra = PETSC_FALSE; 14590598e1a2SLisandro Dalcin PetscMPIInt rank; 1460331830f3SMatthew G. Knepley PetscErrorCode ierr; 1461331830f3SMatthew G. Knepley 1462331830f3SMatthew G. Knepley PetscFunctionBegin; 1463ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 1464ef12879bSLisandro Dalcin ierr = PetscObjectOptionsBegin((PetscObject)viewer);CHKERRQ(ierr); 1465ef12879bSLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Gmsh options");CHKERRQ(ierr); 14660598e1a2SLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", hybrid, &hybrid, NULL);CHKERRQ(ierr); 1467ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL);CHKERRQ(ierr); 1468066ea43fSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_highorder","Generate high-order coordinates", "DMPlexCreateGmsh", highOrder, &highOrder, &highOrderSet);CHKERRQ(ierr); 1469066ea43fSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_project", "Project high-order coordinates to a different space", "DMPlexCreateGmsh", project, &project, NULL);CHKERRQ(ierr); 1470ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL);CHKERRQ(ierr); 1471*a45dabc8SMatthew G. Knepley ierr = PetscOptionsBool("-dm_plex_gmsh_use_regions", "Generate labels with region names", "DMPlexCreateGmsh", useregions, &useregions, NULL);CHKERRQ(ierr); 14720598e1a2SLisandro Dalcin ierr = PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL, PETSC_DECIDE);CHKERRQ(ierr); 1473ef12879bSLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 1474ef12879bSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 14750598e1a2SLisandro Dalcin 14760598e1a2SLisandro Dalcin ierr = GmshCellInfoSetUp();CHKERRQ(ierr); 147711c56cb0SLisandro Dalcin 1478331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1479331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 14800598e1a2SLisandro Dalcin ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr); 148111c56cb0SLisandro Dalcin 148211c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 148311c56cb0SLisandro Dalcin 148411c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 14853b46f529SLisandro Dalcin if (binary) { 148611c56cb0SLisandro Dalcin parentviewer = viewer; 148711c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 148811c56cb0SLisandro Dalcin } 148911c56cb0SLisandro Dalcin 1490dd400576SPatrick Sanan if (rank == 0) { 14910598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1492698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1493698a79b9SLisandro Dalcin PetscBool match; 1494331830f3SMatthew G. Knepley 1495580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 1496698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1497698a79b9SLisandro Dalcin gmsh->binary = binary; 1498698a79b9SLisandro Dalcin 14990598e1a2SLisandro Dalcin ierr = GmshMeshCreate(&mesh);CHKERRQ(ierr); 15000598e1a2SLisandro Dalcin 1501698a79b9SLisandro Dalcin /* Read mesh format */ 1502d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1503d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 15040598e1a2SLisandro Dalcin ierr = GmshReadMeshFormat(gmsh);CHKERRQ(ierr); 1505d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr); 150611c56cb0SLisandro Dalcin 1507dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 1508d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1509d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$PhysicalNames", line, &match);CHKERRQ(ierr); 1510dc0ede3bSMatthew G. Knepley if (match) { 15110598e1a2SLisandro Dalcin ierr = GmshExpect(gmsh, "$PhysicalNames", line);CHKERRQ(ierr); 1512*a45dabc8SMatthew G. Knepley ierr = GmshReadPhysicalNames(gmsh, mesh);CHKERRQ(ierr); 1513d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr); 1514698a79b9SLisandro Dalcin /* Initial read for entity section */ 1515d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1516dc0ede3bSMatthew G. Knepley } 151711c56cb0SLisandro Dalcin 1518de78e4feSLisandro Dalcin /* Read entities */ 1519698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 1520d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr); 15210598e1a2SLisandro Dalcin ierr = GmshReadEntities(gmsh, mesh);CHKERRQ(ierr); 1522d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr); 1523698a79b9SLisandro Dalcin /* Initial read for nodes section */ 1524d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1525de78e4feSLisandro Dalcin } 1526de78e4feSLisandro Dalcin 1527698a79b9SLisandro Dalcin /* Read nodes */ 1528d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr); 15290598e1a2SLisandro Dalcin ierr = GmshReadNodes(gmsh, mesh);CHKERRQ(ierr); 1530d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr); 153111c56cb0SLisandro Dalcin 1532698a79b9SLisandro Dalcin /* Read elements */ 1533feb237baSPierre Jolivet ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1534d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr); 15350598e1a2SLisandro Dalcin ierr = GmshReadElements(gmsh, mesh);CHKERRQ(ierr); 1536d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr); 1537de78e4feSLisandro Dalcin 15380598e1a2SLisandro Dalcin /* Read periodic section (OPTIONAL) */ 1539698a79b9SLisandro Dalcin if (periodic) { 1540ef12879bSLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1541ef12879bSLisandro Dalcin ierr = GmshMatch(gmsh, "$Periodic", line, &periodic);CHKERRQ(ierr); 1542ef12879bSLisandro Dalcin } 1543ef12879bSLisandro Dalcin if (periodic) { 1544d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr); 15450598e1a2SLisandro Dalcin ierr = GmshReadPeriodic(gmsh, mesh);CHKERRQ(ierr); 1546d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr); 1547698a79b9SLisandro Dalcin } 1548698a79b9SLisandro Dalcin 1549698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 1550698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 15510598e1a2SLisandro Dalcin ierr = PetscFree(gmsh->nbuf);CHKERRQ(ierr); 15520598e1a2SLisandro Dalcin 15530598e1a2SLisandro Dalcin dim = mesh->dim; 1554066ea43fSLisandro Dalcin order = mesh->order; 15550598e1a2SLisandro Dalcin numNodes = mesh->numNodes; 15560598e1a2SLisandro Dalcin numElems = mesh->numElems; 15570598e1a2SLisandro Dalcin numVerts = mesh->numVerts; 15580598e1a2SLisandro Dalcin numCells = mesh->numCells; 1559066ea43fSLisandro Dalcin 1560066ea43fSLisandro Dalcin { 1561066ea43fSLisandro Dalcin GmshElement *elemA = mesh->numCells > 0 ? mesh->elements : NULL; 1562066ea43fSLisandro Dalcin GmshElement *elemB = elemA ? elemA + mesh->numCells - 1 : NULL; 1563066ea43fSLisandro Dalcin int ptA = elemA ? GmshCellMap[elemA->cellType].polytope : -1; 1564066ea43fSLisandro Dalcin int ptB = elemB ? GmshCellMap[elemB->cellType].polytope : -1; 1565066ea43fSLisandro Dalcin isSimplex = (ptA == GMSH_QUA || ptA == GMSH_HEX) ? PETSC_FALSE : PETSC_TRUE; 1566066ea43fSLisandro Dalcin isHybrid = (ptA == ptB) ? PETSC_FALSE : PETSC_TRUE; 1567066ea43fSLisandro Dalcin hasTetra = (ptA == GMSH_TET) ? PETSC_TRUE : PETSC_FALSE; 1568066ea43fSLisandro Dalcin } 1569698a79b9SLisandro Dalcin } 1570698a79b9SLisandro Dalcin 1571698a79b9SLisandro Dalcin if (parentviewer) { 1572698a79b9SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 1573698a79b9SLisandro Dalcin } 1574698a79b9SLisandro Dalcin 1575066ea43fSLisandro Dalcin { 1576066ea43fSLisandro Dalcin int buf[6]; 1577698a79b9SLisandro Dalcin 1578066ea43fSLisandro Dalcin buf[0] = (int)dim; 1579066ea43fSLisandro Dalcin buf[1] = (int)order; 1580066ea43fSLisandro Dalcin buf[2] = periodic; 1581066ea43fSLisandro Dalcin buf[3] = isSimplex; 1582066ea43fSLisandro Dalcin buf[4] = isHybrid; 1583066ea43fSLisandro Dalcin buf[5] = hasTetra; 1584066ea43fSLisandro Dalcin 1585ffc4695bSBarry Smith ierr = MPI_Bcast(buf, 6, MPI_INT, 0, comm);CHKERRMPI(ierr); 1586066ea43fSLisandro Dalcin 1587066ea43fSLisandro Dalcin dim = buf[0]; 1588066ea43fSLisandro Dalcin order = buf[1]; 1589066ea43fSLisandro Dalcin periodic = buf[2] ? PETSC_TRUE : PETSC_FALSE; 1590066ea43fSLisandro Dalcin isSimplex = buf[3] ? PETSC_TRUE : PETSC_FALSE; 1591066ea43fSLisandro Dalcin isHybrid = buf[4] ? PETSC_TRUE : PETSC_FALSE; 1592066ea43fSLisandro Dalcin hasTetra = buf[5] ? PETSC_TRUE : PETSC_FALSE; 1593dea2a3c7SStefano Zampini } 1594dea2a3c7SStefano Zampini 1595066ea43fSLisandro Dalcin if (!highOrderSet) highOrder = (order > 1) ? PETSC_TRUE : PETSC_FALSE; 1596066ea43fSLisandro Dalcin if (highOrder && isHybrid) SETERRQ(comm, PETSC_ERR_SUP, "No support for discretization on hybrid meshes yet"); 1597066ea43fSLisandro Dalcin 15980598e1a2SLisandro Dalcin /* We do not want this label automatically computed, instead we fill it here */ 15990598e1a2SLisandro Dalcin ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr); 160011c56cb0SLisandro Dalcin 1601a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 16020598e1a2SLisandro Dalcin ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 16030598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 16040598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 16050598e1a2SLisandro Dalcin DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType); 16060598e1a2SLisandro Dalcin if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR; 16070598e1a2SLisandro Dalcin ierr = DMPlexSetConeSize(*dm, cell, elem->numVerts);CHKERRQ(ierr); 16080598e1a2SLisandro Dalcin ierr = DMPlexSetCellType(*dm, cell, ctype);CHKERRQ(ierr); 1609db397164SMichael Lange } 16100598e1a2SLisandro Dalcin for (v = numCells; v < numCells+numVerts; ++v) { 161196ca5757SLisandro Dalcin ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr); 161296ca5757SLisandro Dalcin } 1613331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 161496ca5757SLisandro Dalcin 1615a4bb7517SMichael Lange /* Add cell-vertex connections */ 16160598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 16170598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 16180598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 16190598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 16200598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 16210598e1a2SLisandro Dalcin cone[v] = numCells + vv; 1622db397164SMichael Lange } 16230598e1a2SLisandro Dalcin ierr = DMPlexReorderCell(*dm, cell, cone);CHKERRQ(ierr); 16240598e1a2SLisandro Dalcin ierr = DMPlexSetCone(*dm, cell, cone);CHKERRQ(ierr); 1625a4bb7517SMichael Lange } 162696ca5757SLisandro Dalcin 1627c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1628331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1629331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1630331830f3SMatthew G. Knepley if (interpolate) { 16315fd9971aSMatthew G. Knepley DM idm; 1632331830f3SMatthew G. Knepley 1633331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1634331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 1635331830f3SMatthew G. Knepley *dm = idm; 1636331830f3SMatthew G. Knepley } 1637917266a4SLisandro Dalcin 1638f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1639dd400576SPatrick Sanan if (rank == 0 && usemarker) { 1640d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1641d3f73514SStefano Zampini 1642d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 1643d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1644d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1645d3f73514SStefano Zampini PetscInt suppSize; 1646d3f73514SStefano Zampini 1647d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 1648d3f73514SStefano Zampini if (suppSize == 1) { 1649d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1650d3f73514SStefano Zampini 1651d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1652d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 16537dd454faSLisandro Dalcin ierr = DMSetLabelValue_Fast(*dm, &marker, "marker", cone[p], 1);CHKERRQ(ierr); 1654d3f73514SStefano Zampini } 1655d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1656d3f73514SStefano Zampini } 1657d3f73514SStefano Zampini } 1658d3f73514SStefano Zampini } 165916ad7e67SMichael Lange 1660dd400576SPatrick Sanan if (rank == 0) { 1661*a45dabc8SMatthew G. Knepley const PetscInt Nr = useregions ? mesh->numRegions : 0; 166211c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1663d1a54cefSMatthew G. Knepley 1664*a45dabc8SMatthew G. Knepley ierr = PetscCalloc1(Nr, ®ionSets);CHKERRQ(ierr); 166516ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 16660598e1a2SLisandro Dalcin for (cell = 0, e = 0; e < numElems; ++e) { 16670598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 16686e1dd89aSLawrence Mitchell 16696e1dd89aSLawrence Mitchell /* Create cell sets */ 16700598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) { 16710598e1a2SLisandro Dalcin if (elem->numTags > 0) { 1672*a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1673*a45dabc8SMatthew G. Knepley PetscInt r; 1674*a45dabc8SMatthew G. Knepley 1675*a45dabc8SMatthew G. Knepley ierr = DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", cell, tag);CHKERRQ(ierr); 1676*a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1677*a45dabc8SMatthew G. Knepley if (mesh->regionTags[r] == tag) {ierr = DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], cell, tag);CHKERRQ(ierr);} 1678*a45dabc8SMatthew G. Knepley } 16796e1dd89aSLawrence Mitchell } 1680917266a4SLisandro Dalcin cell++; 16816e1dd89aSLawrence Mitchell } 16820c070f12SSander Arens 16830598e1a2SLisandro Dalcin /* Create face sets */ 16840598e1a2SLisandro Dalcin if (interpolate && elem->dim == dim-1) { 16850598e1a2SLisandro Dalcin PetscInt joinSize; 16860598e1a2SLisandro Dalcin const PetscInt *join = NULL; 1687*a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1688*a45dabc8SMatthew G. Knepley PetscInt r; 1689*a45dabc8SMatthew G. Knepley 16900598e1a2SLisandro Dalcin /* Find the relevant facet with vertex joins */ 16910598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 16920598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 16930598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 16940598e1a2SLisandro Dalcin cone[v] = vStart + vv; 16950598e1a2SLisandro Dalcin } 16960598e1a2SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr); 16970598e1a2SLisandro Dalcin if (joinSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Could not determine Plex facet for Gmsh element %D (Plex cell %D)", elem->id, e); 1698*a45dabc8SMatthew G. Knepley ierr = DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], tag);CHKERRQ(ierr); 1699*a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1700*a45dabc8SMatthew G. Knepley if (mesh->regionTags[r] == tag) {ierr = DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], join[0], tag);CHKERRQ(ierr);} 1701*a45dabc8SMatthew G. Knepley } 17020598e1a2SLisandro Dalcin ierr = DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr); 17030598e1a2SLisandro Dalcin } 17040598e1a2SLisandro Dalcin 17050c070f12SSander Arens /* Create vertex sets */ 17060598e1a2SLisandro Dalcin if (elem->dim == 0) { 17070598e1a2SLisandro Dalcin if (elem->numTags > 0) { 17080598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[0]; 17090598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1710*a45dabc8SMatthew G. Knepley const PetscInt tag = elem->tags[0]; 1711*a45dabc8SMatthew G. Knepley PetscInt r; 1712*a45dabc8SMatthew G. Knepley 1713*a45dabc8SMatthew G. Knepley ierr = DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vStart + vv, tag);CHKERRQ(ierr); 1714*a45dabc8SMatthew G. Knepley for (r = 0; r < Nr; ++r) { 1715*a45dabc8SMatthew G. Knepley if (mesh->regionTags[r] == tag) {ierr = DMSetLabelValue_Fast(*dm, ®ionSets[r], mesh->regionNames[r], vStart + vv, tag);CHKERRQ(ierr);} 17160598e1a2SLisandro Dalcin } 17170598e1a2SLisandro Dalcin } 17180598e1a2SLisandro Dalcin } 17190598e1a2SLisandro Dalcin } 1720*a45dabc8SMatthew G. Knepley ierr = PetscFree(regionSets);CHKERRQ(ierr); 1721*a45dabc8SMatthew G. Knepley } 17220598e1a2SLisandro Dalcin 17237dd454faSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 17247dd454faSLisandro Dalcin enum {n = 4}; 17257dd454faSLisandro Dalcin PetscBool flag[n]; 17267dd454faSLisandro Dalcin 17277dd454faSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 17287dd454faSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 17297dd454faSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 17307dd454faSLisandro Dalcin flag[3] = marker ? PETSC_TRUE : PETSC_FALSE; 17317dd454faSLisandro Dalcin ierr = MPI_Bcast(flag, n, MPIU_BOOL, 0, comm);CHKERRMPI(ierr); 17327dd454faSLisandro Dalcin if (flag[0]) {ierr = DMCreateLabel(*dm, "Cell Sets");CHKERRQ(ierr);} 17337dd454faSLisandro Dalcin if (flag[1]) {ierr = DMCreateLabel(*dm, "Face Sets");CHKERRQ(ierr);} 17347dd454faSLisandro Dalcin if (flag[2]) {ierr = DMCreateLabel(*dm, "Vertex Sets");CHKERRQ(ierr);} 17357dd454faSLisandro Dalcin if (flag[3]) {ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr);} 17367dd454faSLisandro Dalcin } 17377dd454faSLisandro Dalcin 17380598e1a2SLisandro Dalcin if (periodic) { 17390598e1a2SLisandro Dalcin ierr = PetscBTCreate(numVerts, &periodicVerts);CHKERRQ(ierr); 17400598e1a2SLisandro Dalcin for (n = 0; n < numNodes; ++n) { 17410598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) { 17420598e1a2SLisandro Dalcin if (PetscUnlikely(mesh->periodMap[n] != n)) { 17430598e1a2SLisandro Dalcin PetscInt m = mesh->periodMap[n]; 17440598e1a2SLisandro Dalcin ierr = PetscBTSet(periodicVerts, mesh->vertexMap[n]);CHKERRQ(ierr); 17450598e1a2SLisandro Dalcin ierr = PetscBTSet(periodicVerts, mesh->vertexMap[m]);CHKERRQ(ierr); 17460598e1a2SLisandro Dalcin } 17470598e1a2SLisandro Dalcin } 17480598e1a2SLisandro Dalcin } 17490598e1a2SLisandro Dalcin ierr = PetscBTCreate(numCells, &periodicCells);CHKERRQ(ierr); 17500598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 17510598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 17520598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 17530598e1a2SLisandro Dalcin PetscInt nn = elem->nodes[v]; 17540598e1a2SLisandro Dalcin PetscInt vv = mesh->vertexMap[nn]; 17550598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) { 17560598e1a2SLisandro Dalcin ierr = PetscBTSet(periodicCells, cell);CHKERRQ(ierr); break; 17570c070f12SSander Arens } 17580c070f12SSander Arens } 17590c070f12SSander Arens } 176016ad7e67SMichael Lange } 176116ad7e67SMichael Lange 1762066ea43fSLisandro Dalcin /* Setup coordinate DM */ 17630598e1a2SLisandro Dalcin if (coordDim < 0) coordDim = dim; 17640598e1a2SLisandro Dalcin ierr = DMSetCoordinateDim(*dm, coordDim);CHKERRQ(ierr); 1765066ea43fSLisandro Dalcin ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 1766066ea43fSLisandro Dalcin if (highOrder) { 1767066ea43fSLisandro Dalcin PetscFE fe; 1768066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1769066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_EQUISPACED; 1770066ea43fSLisandro Dalcin 1771066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1772066ea43fSLisandro Dalcin 1773066ea43fSLisandro Dalcin ierr = GmshCreateFE(comm, NULL, isSimplex, continuity, nodeType, dim, coordDim, order, &fe);CHKERRQ(ierr); 1774066ea43fSLisandro Dalcin ierr = PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_fe_view");CHKERRQ(ierr); 1775066ea43fSLisandro Dalcin ierr = DMSetField(cdm, 0, NULL, (PetscObject) fe);CHKERRQ(ierr); 1776066ea43fSLisandro Dalcin ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 1777066ea43fSLisandro Dalcin ierr = DMCreateDS(cdm);CHKERRQ(ierr); 1778066ea43fSLisandro Dalcin } 1779066ea43fSLisandro Dalcin 1780066ea43fSLisandro Dalcin /* Create coordinates */ 1781066ea43fSLisandro Dalcin if (highOrder) { 1782066ea43fSLisandro Dalcin 1783066ea43fSLisandro Dalcin PetscInt maxDof = GmshNumNodes_HEX(order)*coordDim; 1784066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1785066ea43fSLisandro Dalcin PetscSection section; 1786066ea43fSLisandro Dalcin PetscScalar *cellCoords; 1787066ea43fSLisandro Dalcin 1788066ea43fSLisandro Dalcin ierr = DMSetLocalSection(cdm, NULL);CHKERRQ(ierr); 1789066ea43fSLisandro Dalcin ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 1790066ea43fSLisandro Dalcin ierr = PetscSectionClone(coordSection, §ion);CHKERRQ(ierr); 1791066ea43fSLisandro Dalcin ierr = DMPlexSetClosurePermutationTensor(cdm, 0, section);CHKERRQ(ierr); /* XXX Implement DMPlexSetClosurePermutationLexicographic() */ 1792066ea43fSLisandro Dalcin 1793066ea43fSLisandro Dalcin ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 1794066ea43fSLisandro Dalcin ierr = PetscMalloc1(maxDof, &cellCoords);CHKERRQ(ierr); 1795066ea43fSLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1796066ea43fSLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1797066ea43fSLisandro Dalcin const int *lexorder = GmshCellMap[elem->cellType].lexorder(); 1798066ea43fSLisandro Dalcin for (n = 0; n < elem->numNodes; ++n) { 1799066ea43fSLisandro Dalcin const PetscInt node = elem->nodes[lexorder[n]]; 1800066ea43fSLisandro Dalcin for (d = 0; d < coordDim; ++d) 1801066ea43fSLisandro Dalcin cellCoords[n*coordDim+d] = (PetscReal) coords[node*3+d]; 1802066ea43fSLisandro Dalcin } 1803066ea43fSLisandro Dalcin ierr = DMPlexVecSetClosure(cdm, section, coordinates, cell, cellCoords, INSERT_VALUES);CHKERRQ(ierr); 1804066ea43fSLisandro Dalcin } 1805066ea43fSLisandro Dalcin ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 1806066ea43fSLisandro Dalcin ierr = PetscFree(cellCoords);CHKERRQ(ierr); 1807066ea43fSLisandro Dalcin 1808066ea43fSLisandro Dalcin } else { 1809066ea43fSLisandro Dalcin 1810066ea43fSLisandro Dalcin PetscInt *nodeMap; 1811066ea43fSLisandro Dalcin double *coords = mesh ? mesh->nodelist->xyz : NULL; 1812066ea43fSLisandro Dalcin PetscScalar *pointCoords; 1813066ea43fSLisandro Dalcin 1814066ea43fSLisandro Dalcin ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 1815331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 18160598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldComponents(coordSection, 0, coordDim);CHKERRQ(ierr); 18170598e1a2SLisandro Dalcin if (periodic) { /* we need to localize coordinates on cells */ 18180598e1a2SLisandro Dalcin ierr = PetscSectionSetChart(coordSection, 0, numCells+numVerts);CHKERRQ(ierr); 1819f45c9589SStefano Zampini } else { 18200598e1a2SLisandro Dalcin ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVerts);CHKERRQ(ierr); 1821f45c9589SStefano Zampini } 18220598e1a2SLisandro Dalcin if (periodic) { 18230598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 18240598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 18250598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 18260598e1a2SLisandro Dalcin PetscInt dof = elem->numVerts * coordDim; 18270598e1a2SLisandro Dalcin ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr); 18280598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr); 1829f45c9589SStefano Zampini } 1830f45c9589SStefano Zampini } 1831f45c9589SStefano Zampini } 18320598e1a2SLisandro Dalcin for (v = numCells; v < numCells+numVerts; ++v) { 18330598e1a2SLisandro Dalcin ierr = PetscSectionSetDof(coordSection, v, coordDim);CHKERRQ(ierr); 18340598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, v, 0, coordDim);CHKERRQ(ierr); 18350598e1a2SLisandro Dalcin } 1836331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1837066ea43fSLisandro Dalcin 1838066ea43fSLisandro Dalcin ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 1839066ea43fSLisandro Dalcin ierr = VecGetArray(coordinates, &pointCoords);CHKERRQ(ierr); 18400598e1a2SLisandro Dalcin if (periodic) { 18410598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 18420598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 18430598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 18440598e1a2SLisandro Dalcin PetscInt off, node; 18450598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 18460598e1a2SLisandro Dalcin cone[v] = elem->nodes[v]; 1847066ea43fSLisandro Dalcin ierr = DMPlexReorderCell(cdm, cell, cone);CHKERRQ(ierr); 18480598e1a2SLisandro Dalcin ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr); 18490598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 18500598e1a2SLisandro Dalcin for (node = cone[v], d = 0; d < coordDim; ++d) 1851066ea43fSLisandro Dalcin pointCoords[off++] = (PetscReal) coords[node*3+d]; 1852331830f3SMatthew G. Knepley } 1853331830f3SMatthew G. Knepley } 1854331830f3SMatthew G. Knepley } 1855066ea43fSLisandro Dalcin ierr = PetscMalloc1(numVerts, &nodeMap);CHKERRQ(ierr); 18560598e1a2SLisandro Dalcin for (n = 0; n < numNodes; n++) 18570598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) 1858066ea43fSLisandro Dalcin nodeMap[mesh->vertexMap[n]] = n; 18590598e1a2SLisandro Dalcin for (v = 0; v < numVerts; ++v) { 1860066ea43fSLisandro Dalcin PetscInt off, node = nodeMap[v]; 18610598e1a2SLisandro Dalcin ierr = PetscSectionGetOffset(coordSection, numCells + v, &off);CHKERRQ(ierr); 18620598e1a2SLisandro Dalcin for (d = 0; d < coordDim; ++d) 1863066ea43fSLisandro Dalcin pointCoords[off+d] = (PetscReal) coords[node*3+d]; 18640598e1a2SLisandro Dalcin } 1865066ea43fSLisandro Dalcin ierr = PetscFree(nodeMap);CHKERRQ(ierr); 1866066ea43fSLisandro Dalcin ierr = VecRestoreArray(coordinates, &pointCoords);CHKERRQ(ierr); 1867066ea43fSLisandro Dalcin 1868066ea43fSLisandro Dalcin } 1869066ea43fSLisandro Dalcin 1870066ea43fSLisandro Dalcin ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1871066ea43fSLisandro Dalcin ierr = VecSetBlockSize(coordinates, coordDim);CHKERRQ(ierr); 1872331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 18730598e1a2SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 187411c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 187511c56cb0SLisandro Dalcin 18760598e1a2SLisandro Dalcin ierr = GmshMeshDestroy(&mesh);CHKERRQ(ierr); 18770598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&periodicVerts);CHKERRQ(ierr); 18780598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&periodicCells);CHKERRQ(ierr); 187911c56cb0SLisandro Dalcin 1880066ea43fSLisandro Dalcin if (highOrder && project) { 1881066ea43fSLisandro Dalcin PetscFE fe; 1882066ea43fSLisandro Dalcin const char prefix[] = "dm_plex_gmsh_project_"; 1883066ea43fSLisandro Dalcin PetscBool continuity = periodic ? PETSC_FALSE : PETSC_TRUE; 1884066ea43fSLisandro Dalcin PetscDTNodeType nodeType = PETSCDTNODES_GAUSSJACOBI; 1885066ea43fSLisandro Dalcin 1886066ea43fSLisandro Dalcin if (isSimplex) continuity = PETSC_FALSE; /* XXX FIXME Requires DMPlexSetClosurePermutationLexicographic() */ 1887066ea43fSLisandro Dalcin 1888066ea43fSLisandro Dalcin ierr = GmshCreateFE(comm, prefix, isSimplex, continuity, nodeType, dim, coordDim, order, &fe);CHKERRQ(ierr); 1889066ea43fSLisandro Dalcin ierr = PetscFEViewFromOptions(fe, NULL, "-dm_plex_gmsh_project_fe_view");CHKERRQ(ierr); 1890066ea43fSLisandro Dalcin ierr = DMProjectCoordinates(*dm, fe);CHKERRQ(ierr); 1891066ea43fSLisandro Dalcin ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 1892066ea43fSLisandro Dalcin } 1893066ea43fSLisandro Dalcin 18940598e1a2SLisandro Dalcin ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr); 1895331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1896331830f3SMatthew G. Knepley } 1897