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 5de78e4feSLisandro Dalcin typedef struct { 6*0598e1a2SLisandro Dalcin int cellType; 7*0598e1a2SLisandro Dalcin DMPolytopeType polytope; 8*0598e1a2SLisandro Dalcin int dim; 9*0598e1a2SLisandro Dalcin int numVerts; 10*0598e1a2SLisandro Dalcin int order; 11*0598e1a2SLisandro Dalcin int numNodes; 12*0598e1a2SLisandro Dalcin } GmshCellInfo; 13*0598e1a2SLisandro Dalcin 14*0598e1a2SLisandro Dalcin #define DM_POLYTOPE_VERTEX DM_POLYTOPE_POINT 15*0598e1a2SLisandro Dalcin #define DM_POLYTOPE_PYRAMID DM_POLYTOPE_UNKNOWN 16*0598e1a2SLisandro Dalcin 17*0598e1a2SLisandro Dalcin static const GmshCellInfo GmshCellTable[] = { 18*0598e1a2SLisandro Dalcin { 15, DM_POLYTOPE_VERTEX, 0, 1, 0, 1}, 19*0598e1a2SLisandro Dalcin 20*0598e1a2SLisandro Dalcin { 1, DM_POLYTOPE_SEGMENT, 1, 2, 1, 2}, 21*0598e1a2SLisandro Dalcin { 8, DM_POLYTOPE_SEGMENT, 1, 2, 2, 3}, 22*0598e1a2SLisandro Dalcin { 26, DM_POLYTOPE_SEGMENT, 1, 2, 3, 4}, 23*0598e1a2SLisandro Dalcin { 27, DM_POLYTOPE_SEGMENT, 1, 2, 4, 5}, 24*0598e1a2SLisandro Dalcin { 28, DM_POLYTOPE_SEGMENT, 1, 2, 5, 6}, 25*0598e1a2SLisandro Dalcin { 62, DM_POLYTOPE_SEGMENT, 1, 2, 6, 7}, 26*0598e1a2SLisandro Dalcin { 63, DM_POLYTOPE_SEGMENT, 1, 2, 7, 8}, 27*0598e1a2SLisandro Dalcin { 64, DM_POLYTOPE_SEGMENT, 1, 2, 8, 9}, 28*0598e1a2SLisandro Dalcin { 65, DM_POLYTOPE_SEGMENT, 1, 2, 9, 10}, 29*0598e1a2SLisandro Dalcin { 66, DM_POLYTOPE_SEGMENT, 1, 2, 10, 11}, 30*0598e1a2SLisandro Dalcin 31*0598e1a2SLisandro Dalcin { 2, DM_POLYTOPE_TRIANGLE, 2, 3, 1, 3}, 32*0598e1a2SLisandro Dalcin { 9, DM_POLYTOPE_TRIANGLE, 2, 3, 2, 6}, 33*0598e1a2SLisandro Dalcin { 21, DM_POLYTOPE_TRIANGLE, 2, 3, 3, 10}, 34*0598e1a2SLisandro Dalcin { 23, DM_POLYTOPE_TRIANGLE, 2, 3, 4, 15}, 35*0598e1a2SLisandro Dalcin { 25, DM_POLYTOPE_TRIANGLE, 2, 3, 5, 21}, 36*0598e1a2SLisandro Dalcin { 42, DM_POLYTOPE_TRIANGLE, 2, 3, 6, 28}, 37*0598e1a2SLisandro Dalcin { 43, DM_POLYTOPE_TRIANGLE, 2, 3, 7, 36}, 38*0598e1a2SLisandro Dalcin { 44, DM_POLYTOPE_TRIANGLE, 2, 3, 8, 45}, 39*0598e1a2SLisandro Dalcin { 45, DM_POLYTOPE_TRIANGLE, 2, 3, 9, 55}, 40*0598e1a2SLisandro Dalcin { 46, DM_POLYTOPE_TRIANGLE, 2, 3, 10, 66}, 41*0598e1a2SLisandro Dalcin 42*0598e1a2SLisandro Dalcin { 3, DM_POLYTOPE_QUADRILATERAL, 2, 4, 1, 4}, 43*0598e1a2SLisandro Dalcin { 10, DM_POLYTOPE_QUADRILATERAL, 2, 4, 2, 9}, 44*0598e1a2SLisandro Dalcin { 36, DM_POLYTOPE_QUADRILATERAL, 2, 4, 3, 16}, 45*0598e1a2SLisandro Dalcin { 37, DM_POLYTOPE_QUADRILATERAL, 2, 4, 4, 25}, 46*0598e1a2SLisandro Dalcin { 38, DM_POLYTOPE_QUADRILATERAL, 2, 4, 5, 36}, 47*0598e1a2SLisandro Dalcin { 47, DM_POLYTOPE_QUADRILATERAL, 2, 4, 6, 49}, 48*0598e1a2SLisandro Dalcin { 48, DM_POLYTOPE_QUADRILATERAL, 2, 4, 7, 64}, 49*0598e1a2SLisandro Dalcin { 49, DM_POLYTOPE_QUADRILATERAL, 2, 4, 8, 81}, 50*0598e1a2SLisandro Dalcin { 50, DM_POLYTOPE_QUADRILATERAL, 2, 4, 9, 100}, 51*0598e1a2SLisandro Dalcin { 51, DM_POLYTOPE_QUADRILATERAL, 2, 4, 10, 121}, 52*0598e1a2SLisandro Dalcin 53*0598e1a2SLisandro Dalcin { 4, DM_POLYTOPE_TETRAHEDRON, 3, 4, 1, 4}, 54*0598e1a2SLisandro Dalcin { 11, DM_POLYTOPE_TETRAHEDRON, 3, 4, 2, 10}, 55*0598e1a2SLisandro Dalcin { 29, DM_POLYTOPE_TETRAHEDRON, 3, 4, 3, 20}, 56*0598e1a2SLisandro Dalcin { 30, DM_POLYTOPE_TETRAHEDRON, 3, 4, 4, 35}, 57*0598e1a2SLisandro Dalcin { 31, DM_POLYTOPE_TETRAHEDRON, 3, 4, 5, 56}, 58*0598e1a2SLisandro Dalcin { 71, DM_POLYTOPE_TETRAHEDRON, 3, 4, 6, 84}, 59*0598e1a2SLisandro Dalcin { 72, DM_POLYTOPE_TETRAHEDRON, 3, 4, 7, 120}, 60*0598e1a2SLisandro Dalcin { 73, DM_POLYTOPE_TETRAHEDRON, 3, 4, 8, 165}, 61*0598e1a2SLisandro Dalcin { 74, DM_POLYTOPE_TETRAHEDRON, 3, 4, 9, 220}, 62*0598e1a2SLisandro Dalcin { 75, DM_POLYTOPE_TETRAHEDRON, 3, 4, 10, 286}, 63*0598e1a2SLisandro Dalcin 64*0598e1a2SLisandro Dalcin { 5, DM_POLYTOPE_HEXAHEDRON, 3, 8, 1, 8}, 65*0598e1a2SLisandro Dalcin { 12, DM_POLYTOPE_HEXAHEDRON, 3, 8, 2, 27}, 66*0598e1a2SLisandro Dalcin { 92, DM_POLYTOPE_HEXAHEDRON, 3, 8, 3, 64}, 67*0598e1a2SLisandro Dalcin { 93, DM_POLYTOPE_HEXAHEDRON, 3, 8, 4, 125}, 68*0598e1a2SLisandro Dalcin { 94, DM_POLYTOPE_HEXAHEDRON, 3, 8, 5, 216}, 69*0598e1a2SLisandro Dalcin { 95, DM_POLYTOPE_HEXAHEDRON, 3, 8, 6, 343}, 70*0598e1a2SLisandro Dalcin { 96, DM_POLYTOPE_HEXAHEDRON, 3, 8, 7, 512}, 71*0598e1a2SLisandro Dalcin { 97, DM_POLYTOPE_HEXAHEDRON, 3, 8, 8, 729}, 72*0598e1a2SLisandro Dalcin { 98, DM_POLYTOPE_HEXAHEDRON, 3, 8, 9, 1000}, 73*0598e1a2SLisandro Dalcin 74*0598e1a2SLisandro Dalcin { 6, DM_POLYTOPE_TRI_PRISM, 3, 6, 1, 6}, 75*0598e1a2SLisandro Dalcin { 13, DM_POLYTOPE_TRI_PRISM, 3, 6, 2, 18}, 76*0598e1a2SLisandro Dalcin { 90, DM_POLYTOPE_TRI_PRISM, 3, 6, 3, 40}, 77*0598e1a2SLisandro Dalcin { 91, DM_POLYTOPE_TRI_PRISM, 3, 6, 4, 75}, 78*0598e1a2SLisandro Dalcin {106, DM_POLYTOPE_TRI_PRISM, 3, 6, 5, 126}, 79*0598e1a2SLisandro Dalcin {107, DM_POLYTOPE_TRI_PRISM, 3, 6, 6, 196}, 80*0598e1a2SLisandro Dalcin {108, DM_POLYTOPE_TRI_PRISM, 3, 6, 7, 288}, 81*0598e1a2SLisandro Dalcin {109, DM_POLYTOPE_TRI_PRISM, 3, 6, 8, 405}, 82*0598e1a2SLisandro Dalcin {110, DM_POLYTOPE_TRI_PRISM, 3, 6, 9, 550}, 83*0598e1a2SLisandro Dalcin 84*0598e1a2SLisandro Dalcin { 7, DM_POLYTOPE_PYRAMID, 3, 5, 1, 5}, 85*0598e1a2SLisandro Dalcin { 14, DM_POLYTOPE_PYRAMID, 3, 5, 2, 14}, 86*0598e1a2SLisandro Dalcin {118, DM_POLYTOPE_PYRAMID, 3, 5, 3, 30}, 87*0598e1a2SLisandro Dalcin {119, DM_POLYTOPE_PYRAMID, 3, 5, 4, 55}, 88*0598e1a2SLisandro Dalcin {120, DM_POLYTOPE_PYRAMID, 3, 5, 5, 91}, 89*0598e1a2SLisandro Dalcin {121, DM_POLYTOPE_PYRAMID, 3, 5, 6, 140}, 90*0598e1a2SLisandro Dalcin {122, DM_POLYTOPE_PYRAMID, 3, 5, 7, 204}, 91*0598e1a2SLisandro Dalcin {123, DM_POLYTOPE_PYRAMID, 3, 5, 8, 285}, 92*0598e1a2SLisandro Dalcin {124, DM_POLYTOPE_PYRAMID, 3, 5, 9, 385}, 93*0598e1a2SLisandro Dalcin 94*0598e1a2SLisandro Dalcin #if 0 95*0598e1a2SLisandro Dalcin { 20, DM_POLYTOPE_TRIANGLE, 2, 3, 3, 9}, 96*0598e1a2SLisandro Dalcin { 16, DM_POLYTOPE_QUADRILATERAL, 2, 4, 2, 8}, 97*0598e1a2SLisandro Dalcin { 17, DM_POLYTOPE_HEXAHEDRON, 3, 8, 2, 20}, 98*0598e1a2SLisandro Dalcin { 18, DM_POLYTOPE_TRI_PRISM, 3, 6, 2, 15}, 99*0598e1a2SLisandro Dalcin { 19, DM_POLYTOPE_PYRAMID, 3, 5, 2, 13}, 100*0598e1a2SLisandro Dalcin #endif 101*0598e1a2SLisandro Dalcin }; 102*0598e1a2SLisandro Dalcin 103*0598e1a2SLisandro Dalcin static GmshCellInfo GmshCellMap[150]; 104*0598e1a2SLisandro Dalcin 105*0598e1a2SLisandro Dalcin static PetscErrorCode GmshCellInfoSetUp(void) 106*0598e1a2SLisandro Dalcin { 107*0598e1a2SLisandro Dalcin size_t i, n; 108*0598e1a2SLisandro Dalcin static PetscBool called = PETSC_FALSE; 109*0598e1a2SLisandro Dalcin 110*0598e1a2SLisandro Dalcin if (called) return 0; 111*0598e1a2SLisandro Dalcin PetscFunctionBegin; 112*0598e1a2SLisandro Dalcin called = PETSC_TRUE; 113*0598e1a2SLisandro Dalcin n = sizeof(GmshCellMap)/sizeof(GmshCellMap[0]); 114*0598e1a2SLisandro Dalcin for (i = 0; i < n; ++i) { 115*0598e1a2SLisandro Dalcin GmshCellMap[i].cellType = -1; 116*0598e1a2SLisandro Dalcin GmshCellMap[i].polytope = DM_POLYTOPE_UNKNOWN; 117*0598e1a2SLisandro Dalcin } 118*0598e1a2SLisandro Dalcin n = sizeof(GmshCellTable)/sizeof(GmshCellTable[0]); 119*0598e1a2SLisandro Dalcin for (i = 0; i < n; ++i) GmshCellMap[GmshCellTable[i].cellType] = GmshCellTable[i]; 120*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 121*0598e1a2SLisandro Dalcin } 122*0598e1a2SLisandro Dalcin 123*0598e1a2SLisandro Dalcin #define GmshCellTypeCheck(ct) 0; do { \ 124*0598e1a2SLisandro Dalcin const int _ct_ = (int)ct; \ 125*0598e1a2SLisandro Dalcin if (_ct_ < 0 || _ct_ >= (int)(sizeof(GmshCellMap)/sizeof(GmshCellMap[0]))) \ 126*0598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid Gmsh element type %d", _ct_); \ 127*0598e1a2SLisandro Dalcin if (GmshCellMap[_ct_].cellType != _ct_) \ 128*0598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 129*0598e1a2SLisandro Dalcin if (GmshCellMap[_ct_].polytope == DM_POLYTOPE_UNKNOWN) \ 130*0598e1a2SLisandro Dalcin SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported Gmsh element type %d", _ct_); \ 131*0598e1a2SLisandro Dalcin } while(0) 132*0598e1a2SLisandro Dalcin 133*0598e1a2SLisandro Dalcin 134*0598e1a2SLisandro Dalcin typedef struct { 135698a79b9SLisandro Dalcin PetscViewer viewer; 136698a79b9SLisandro Dalcin int fileFormat; 137698a79b9SLisandro Dalcin int dataSize; 138698a79b9SLisandro Dalcin PetscBool binary; 139698a79b9SLisandro Dalcin PetscBool byteSwap; 140698a79b9SLisandro Dalcin size_t wlen; 141698a79b9SLisandro Dalcin void *wbuf; 142698a79b9SLisandro Dalcin size_t slen; 143698a79b9SLisandro Dalcin void *sbuf; 144*0598e1a2SLisandro Dalcin PetscInt *nbuf; 145*0598e1a2SLisandro Dalcin PetscInt nodeStart; 146*0598e1a2SLisandro Dalcin PetscInt nodeEnd; 14733a088b6SMatthew G. Knepley PetscInt *nodeMap; 148698a79b9SLisandro Dalcin } GmshFile; 149de78e4feSLisandro Dalcin 150698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferGet(GmshFile *gmsh, size_t count, size_t eltsize, void *buf) 151de78e4feSLisandro Dalcin { 152de78e4feSLisandro Dalcin size_t size = count * eltsize; 15311c56cb0SLisandro Dalcin PetscErrorCode ierr; 15411c56cb0SLisandro Dalcin 15511c56cb0SLisandro Dalcin PetscFunctionBegin; 156698a79b9SLisandro Dalcin if (gmsh->wlen < size) { 157698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 158698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->wbuf);CHKERRQ(ierr); 159698a79b9SLisandro Dalcin gmsh->wlen = size; 160de78e4feSLisandro Dalcin } 161698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->wbuf : NULL; 162de78e4feSLisandro Dalcin PetscFunctionReturn(0); 163de78e4feSLisandro Dalcin } 164de78e4feSLisandro Dalcin 165698a79b9SLisandro Dalcin static PetscErrorCode GmshBufferSizeGet(GmshFile *gmsh, size_t count, void *buf) 166698a79b9SLisandro Dalcin { 167698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 168698a79b9SLisandro Dalcin size_t size = count * dataSize; 169698a79b9SLisandro Dalcin PetscErrorCode ierr; 170698a79b9SLisandro Dalcin 171698a79b9SLisandro Dalcin PetscFunctionBegin; 172698a79b9SLisandro Dalcin if (gmsh->slen < size) { 173698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 174698a79b9SLisandro Dalcin ierr = PetscMalloc(size, &gmsh->sbuf);CHKERRQ(ierr); 175698a79b9SLisandro Dalcin gmsh->slen = size; 176698a79b9SLisandro Dalcin } 177698a79b9SLisandro Dalcin *(void**)buf = size ? gmsh->sbuf : NULL; 178698a79b9SLisandro Dalcin PetscFunctionReturn(0); 179698a79b9SLisandro Dalcin } 180698a79b9SLisandro Dalcin 181698a79b9SLisandro Dalcin static PetscErrorCode GmshRead(GmshFile *gmsh, void *buf, PetscInt count, PetscDataType dtype) 182de78e4feSLisandro Dalcin { 183de78e4feSLisandro Dalcin PetscErrorCode ierr; 184de78e4feSLisandro Dalcin PetscFunctionBegin; 185698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, dtype);CHKERRQ(ierr); 186698a79b9SLisandro Dalcin if (gmsh->byteSwap) {ierr = PetscByteSwap(buf, dtype, count);CHKERRQ(ierr);} 187698a79b9SLisandro Dalcin PetscFunctionReturn(0); 188698a79b9SLisandro Dalcin } 189698a79b9SLisandro Dalcin 190698a79b9SLisandro Dalcin static PetscErrorCode GmshReadString(GmshFile *gmsh, char *buf, PetscInt count) 191698a79b9SLisandro Dalcin { 192698a79b9SLisandro Dalcin PetscErrorCode ierr; 193698a79b9SLisandro Dalcin PetscFunctionBegin; 194698a79b9SLisandro Dalcin ierr = PetscViewerRead(gmsh->viewer, buf, count, NULL, PETSC_STRING);CHKERRQ(ierr); 195698a79b9SLisandro Dalcin PetscFunctionReturn(0); 196698a79b9SLisandro Dalcin } 197698a79b9SLisandro Dalcin 198d6689ca9SLisandro Dalcin static PetscErrorCode GmshMatch(PETSC_UNUSED GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN], PetscBool *match) 199d6689ca9SLisandro Dalcin { 200d6689ca9SLisandro Dalcin PetscErrorCode ierr; 201d6689ca9SLisandro Dalcin PetscFunctionBegin; 202d6689ca9SLisandro Dalcin ierr = PetscStrcmp(line, Section, match);CHKERRQ(ierr); 203d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 204d6689ca9SLisandro Dalcin } 205d6689ca9SLisandro Dalcin 206d6689ca9SLisandro Dalcin static PetscErrorCode GmshExpect(GmshFile *gmsh, const char Section[], char line[PETSC_MAX_PATH_LEN]) 207d6689ca9SLisandro Dalcin { 208d6689ca9SLisandro Dalcin PetscBool match; 209d6689ca9SLisandro Dalcin PetscErrorCode ierr; 210d6689ca9SLisandro Dalcin 211d6689ca9SLisandro Dalcin PetscFunctionBegin; 212d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, Section, line, &match);CHKERRQ(ierr); 213*0598e1a2SLisandro Dalcin if (!match) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file, expecting %s",Section); 214d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 215d6689ca9SLisandro Dalcin } 216d6689ca9SLisandro Dalcin 217d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadSection(GmshFile *gmsh, char line[PETSC_MAX_PATH_LEN]) 218d6689ca9SLisandro Dalcin { 219d6689ca9SLisandro Dalcin PetscBool match; 220d6689ca9SLisandro Dalcin PetscErrorCode ierr; 221d6689ca9SLisandro Dalcin 222d6689ca9SLisandro Dalcin PetscFunctionBegin; 223d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 224d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 225d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$Comments", line, &match);CHKERRQ(ierr); 226d6689ca9SLisandro Dalcin if (!match) break; 227d6689ca9SLisandro Dalcin while (PETSC_TRUE) { 228d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 229d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$EndComments", line, &match);CHKERRQ(ierr); 230d6689ca9SLisandro Dalcin if (match) break; 231d6689ca9SLisandro Dalcin } 232d6689ca9SLisandro Dalcin } 233d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 234d6689ca9SLisandro Dalcin } 235d6689ca9SLisandro Dalcin 236d6689ca9SLisandro Dalcin static PetscErrorCode GmshReadEndSection(GmshFile *gmsh, const char EndSection[], char line[PETSC_MAX_PATH_LEN]) 237d6689ca9SLisandro Dalcin { 238d6689ca9SLisandro Dalcin PetscErrorCode ierr; 239d6689ca9SLisandro Dalcin PetscFunctionBegin; 240d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 241d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, EndSection, line);CHKERRQ(ierr); 242d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 243d6689ca9SLisandro Dalcin } 244d6689ca9SLisandro Dalcin 245698a79b9SLisandro Dalcin static PetscErrorCode GmshReadSize(GmshFile *gmsh, PetscInt *buf, PetscInt count) 246698a79b9SLisandro Dalcin { 247698a79b9SLisandro Dalcin PetscInt i; 248698a79b9SLisandro Dalcin size_t dataSize = (size_t)gmsh->dataSize; 249698a79b9SLisandro Dalcin PetscErrorCode ierr; 250698a79b9SLisandro Dalcin 251698a79b9SLisandro Dalcin PetscFunctionBegin; 252698a79b9SLisandro Dalcin if (dataSize == sizeof(PetscInt)) { 253698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_INT);CHKERRQ(ierr); 254698a79b9SLisandro Dalcin } else if (dataSize == sizeof(int)) { 255698a79b9SLisandro Dalcin int *ibuf = NULL; 25689d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 257698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 258698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 259698a79b9SLisandro Dalcin } else if (dataSize == sizeof(long)) { 260698a79b9SLisandro Dalcin long *ibuf = NULL; 26189d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 262698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_LONG);CHKERRQ(ierr); 263698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 264698a79b9SLisandro Dalcin } else if (dataSize == sizeof(PetscInt64)) { 265698a79b9SLisandro Dalcin PetscInt64 *ibuf = NULL; 26689d5b078SSatish Balay ierr = GmshBufferSizeGet(gmsh, count, &ibuf);CHKERRQ(ierr); 267698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, ibuf, count, PETSC_INT64);CHKERRQ(ierr); 268698a79b9SLisandro Dalcin for (i = 0; i < count; ++i) buf[i] = (PetscInt)ibuf[i]; 269698a79b9SLisandro Dalcin } 270698a79b9SLisandro Dalcin PetscFunctionReturn(0); 271698a79b9SLisandro Dalcin } 272698a79b9SLisandro Dalcin 273698a79b9SLisandro Dalcin static PetscErrorCode GmshReadInt(GmshFile *gmsh, int *buf, PetscInt count) 274698a79b9SLisandro Dalcin { 275698a79b9SLisandro Dalcin PetscErrorCode ierr; 276698a79b9SLisandro Dalcin PetscFunctionBegin; 277698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_ENUM);CHKERRQ(ierr); 278698a79b9SLisandro Dalcin PetscFunctionReturn(0); 279698a79b9SLisandro Dalcin } 280698a79b9SLisandro Dalcin 281698a79b9SLisandro Dalcin static PetscErrorCode GmshReadDouble(GmshFile *gmsh, double *buf, PetscInt count) 282698a79b9SLisandro Dalcin { 283698a79b9SLisandro Dalcin PetscErrorCode ierr; 284698a79b9SLisandro Dalcin PetscFunctionBegin; 285698a79b9SLisandro Dalcin ierr = GmshRead(gmsh, buf, count, PETSC_DOUBLE);CHKERRQ(ierr); 286de78e4feSLisandro Dalcin PetscFunctionReturn(0); 287de78e4feSLisandro Dalcin } 288de78e4feSLisandro Dalcin 289de78e4feSLisandro Dalcin typedef struct { 290*0598e1a2SLisandro Dalcin PetscInt id; /* Entity ID */ 291*0598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 292de78e4feSLisandro Dalcin double bbox[6]; /* Bounding box */ 293de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 294de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 295de78e4feSLisandro Dalcin } GmshEntity; 296de78e4feSLisandro Dalcin 297de78e4feSLisandro Dalcin typedef struct { 298730356f6SLisandro Dalcin GmshEntity *entity[4]; 299730356f6SLisandro Dalcin PetscHMapI entityMap[4]; 300730356f6SLisandro Dalcin } GmshEntities; 301730356f6SLisandro Dalcin 302698a79b9SLisandro Dalcin static PetscErrorCode GmshEntitiesCreate(PetscInt count[4], GmshEntities **entities) 303730356f6SLisandro Dalcin { 304730356f6SLisandro Dalcin PetscInt dim; 305730356f6SLisandro Dalcin PetscErrorCode ierr; 306730356f6SLisandro Dalcin 307730356f6SLisandro Dalcin PetscFunctionBegin; 308730356f6SLisandro Dalcin ierr = PetscNew(entities);CHKERRQ(ierr); 309730356f6SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 310730356f6SLisandro Dalcin ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr); 311730356f6SLisandro Dalcin ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 312730356f6SLisandro Dalcin } 313730356f6SLisandro Dalcin PetscFunctionReturn(0); 314730356f6SLisandro Dalcin } 315730356f6SLisandro Dalcin 316*0598e1a2SLisandro Dalcin static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 317*0598e1a2SLisandro Dalcin { 318*0598e1a2SLisandro Dalcin PetscInt dim; 319*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 320*0598e1a2SLisandro Dalcin 321*0598e1a2SLisandro Dalcin PetscFunctionBegin; 322*0598e1a2SLisandro Dalcin if (!*entities) PetscFunctionReturn(0); 323*0598e1a2SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 324*0598e1a2SLisandro Dalcin ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr); 325*0598e1a2SLisandro Dalcin ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 326*0598e1a2SLisandro Dalcin } 327*0598e1a2SLisandro Dalcin ierr = PetscFree((*entities));CHKERRQ(ierr); 328*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 329*0598e1a2SLisandro Dalcin } 330*0598e1a2SLisandro Dalcin 331730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 332730356f6SLisandro Dalcin { 333730356f6SLisandro Dalcin PetscErrorCode ierr; 334*0598e1a2SLisandro Dalcin 335730356f6SLisandro Dalcin PetscFunctionBegin; 336730356f6SLisandro Dalcin ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr); 337730356f6SLisandro Dalcin entities->entity[dim][index].dim = dim; 338730356f6SLisandro Dalcin entities->entity[dim][index].id = eid; 339730356f6SLisandro Dalcin if (entity) *entity = &entities->entity[dim][index]; 340730356f6SLisandro Dalcin PetscFunctionReturn(0); 341730356f6SLisandro Dalcin } 342730356f6SLisandro Dalcin 343730356f6SLisandro Dalcin static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 344730356f6SLisandro Dalcin { 345730356f6SLisandro Dalcin PetscInt index; 346730356f6SLisandro Dalcin PetscErrorCode ierr; 347730356f6SLisandro Dalcin 348730356f6SLisandro Dalcin PetscFunctionBegin; 349730356f6SLisandro Dalcin ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr); 350730356f6SLisandro Dalcin *entity = &entities->entity[dim][index]; 351730356f6SLisandro Dalcin PetscFunctionReturn(0); 352730356f6SLisandro Dalcin } 353730356f6SLisandro Dalcin 354*0598e1a2SLisandro Dalcin typedef struct { 355*0598e1a2SLisandro Dalcin PetscInt *id; /* Node IDs */ 356*0598e1a2SLisandro Dalcin double *xyz; /* Coordinates */ 357*0598e1a2SLisandro Dalcin } GmshNodes; 358*0598e1a2SLisandro Dalcin 359*0598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesCreate(PetscInt count, GmshNodes **nodes) 360730356f6SLisandro Dalcin { 361730356f6SLisandro Dalcin PetscErrorCode ierr; 362730356f6SLisandro Dalcin 363730356f6SLisandro Dalcin PetscFunctionBegin; 364*0598e1a2SLisandro Dalcin ierr = PetscNew(nodes);CHKERRQ(ierr); 365*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(count*1, &(*nodes)->id);CHKERRQ(ierr); 366*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(count*3, &(*nodes)->xyz);CHKERRQ(ierr); 367*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 368730356f6SLisandro Dalcin } 369*0598e1a2SLisandro Dalcin 370*0598e1a2SLisandro Dalcin static PetscErrorCode GmshNodesDestroy(GmshNodes **nodes) 371*0598e1a2SLisandro Dalcin { 372*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 373*0598e1a2SLisandro Dalcin PetscFunctionBegin; 374*0598e1a2SLisandro Dalcin if (!*nodes) PetscFunctionReturn(0); 375*0598e1a2SLisandro Dalcin ierr = PetscFree((*nodes)->id);CHKERRQ(ierr); 376*0598e1a2SLisandro Dalcin ierr = PetscFree((*nodes)->xyz);CHKERRQ(ierr); 377*0598e1a2SLisandro Dalcin ierr = PetscFree((*nodes));CHKERRQ(ierr); 378730356f6SLisandro Dalcin PetscFunctionReturn(0); 379730356f6SLisandro Dalcin } 380730356f6SLisandro Dalcin 381730356f6SLisandro Dalcin typedef struct { 382*0598e1a2SLisandro Dalcin PetscInt id; /* Element ID */ 383*0598e1a2SLisandro Dalcin PetscInt dim; /* Dimension */ 384de78e4feSLisandro Dalcin PetscInt cellType; /* Cell type */ 385*0598e1a2SLisandro Dalcin PetscInt numVerts; /* Size of vertex array */ 386de78e4feSLisandro Dalcin PetscInt numNodes; /* Size of node array */ 387*0598e1a2SLisandro Dalcin PetscInt *nodes; /* Vertex/Node array */ 388de78e4feSLisandro Dalcin PetscInt numTags; /* Size of tag array */ 389de78e4feSLisandro Dalcin int tags[4]; /* Tag array */ 390de78e4feSLisandro Dalcin } GmshElement; 391de78e4feSLisandro Dalcin 392*0598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsCreate(PetscInt count, GmshElement **elements) 393*0598e1a2SLisandro Dalcin { 394*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 395*0598e1a2SLisandro Dalcin 396*0598e1a2SLisandro Dalcin PetscFunctionBegin; 397*0598e1a2SLisandro Dalcin ierr = PetscCalloc1(count, elements);CHKERRQ(ierr); 398*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 399*0598e1a2SLisandro Dalcin } 400*0598e1a2SLisandro Dalcin 401*0598e1a2SLisandro Dalcin static PetscErrorCode GmshElementsDestroy(GmshElement **elements) 402*0598e1a2SLisandro Dalcin { 403*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 404*0598e1a2SLisandro Dalcin 405*0598e1a2SLisandro Dalcin PetscFunctionBegin; 406*0598e1a2SLisandro Dalcin if (!*elements) PetscFunctionReturn(0); 407*0598e1a2SLisandro Dalcin ierr = PetscFree(*elements);CHKERRQ(ierr); 408*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 409*0598e1a2SLisandro Dalcin } 410*0598e1a2SLisandro Dalcin 411*0598e1a2SLisandro Dalcin typedef struct { 412*0598e1a2SLisandro Dalcin PetscInt dim; 413*0598e1a2SLisandro Dalcin GmshEntities *entities; 414*0598e1a2SLisandro Dalcin PetscInt numNodes; 415*0598e1a2SLisandro Dalcin GmshNodes *nodelist; 416*0598e1a2SLisandro Dalcin PetscInt numElems; 417*0598e1a2SLisandro Dalcin GmshElement *elements; 418*0598e1a2SLisandro Dalcin PetscInt numVerts; 419*0598e1a2SLisandro Dalcin PetscInt numCells; 420*0598e1a2SLisandro Dalcin PetscInt *periodMap; 421*0598e1a2SLisandro Dalcin PetscInt *vertexMap; 422*0598e1a2SLisandro Dalcin PetscSegBuffer segbuf; 423*0598e1a2SLisandro Dalcin } GmshMesh; 424*0598e1a2SLisandro Dalcin 425*0598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshCreate(GmshMesh **mesh) 426*0598e1a2SLisandro Dalcin { 427*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 428*0598e1a2SLisandro Dalcin 429*0598e1a2SLisandro Dalcin PetscFunctionBegin; 430*0598e1a2SLisandro Dalcin ierr = PetscNew(mesh);CHKERRQ(ierr); 431*0598e1a2SLisandro Dalcin ierr = PetscSegBufferCreate(sizeof(PetscInt), 0, &(*mesh)->segbuf);CHKERRQ(ierr); 432*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 433*0598e1a2SLisandro Dalcin } 434*0598e1a2SLisandro Dalcin 435*0598e1a2SLisandro Dalcin static PetscErrorCode GmshMeshDestroy(GmshMesh **mesh) 436*0598e1a2SLisandro Dalcin { 437*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 438*0598e1a2SLisandro Dalcin 439*0598e1a2SLisandro Dalcin PetscFunctionBegin; 440*0598e1a2SLisandro Dalcin if (!*mesh) PetscFunctionReturn(0); 441*0598e1a2SLisandro Dalcin ierr = GmshEntitiesDestroy(&(*mesh)->entities);CHKERRQ(ierr); 442*0598e1a2SLisandro Dalcin ierr = GmshNodesDestroy(&(*mesh)->nodelist);CHKERRQ(ierr); 443*0598e1a2SLisandro Dalcin ierr = GmshElementsDestroy(&(*mesh)->elements);CHKERRQ(ierr); 444*0598e1a2SLisandro Dalcin ierr = PetscFree((*mesh)->periodMap);CHKERRQ(ierr); 445*0598e1a2SLisandro Dalcin ierr = PetscFree((*mesh)->vertexMap);CHKERRQ(ierr); 446*0598e1a2SLisandro Dalcin ierr = PetscSegBufferDestroy(&(*mesh)->segbuf);CHKERRQ(ierr); 447*0598e1a2SLisandro Dalcin ierr = PetscFree((*mesh));CHKERRQ(ierr); 448*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 449*0598e1a2SLisandro Dalcin } 450*0598e1a2SLisandro Dalcin 451*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v22(GmshFile *gmsh, GmshMesh *mesh) 452de78e4feSLisandro Dalcin { 453698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 454698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 455de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 456*0598e1a2SLisandro Dalcin int n, num, nid, snum; 457*0598e1a2SLisandro Dalcin GmshNodes *nodes; 458de78e4feSLisandro Dalcin PetscErrorCode ierr; 459de78e4feSLisandro Dalcin 460de78e4feSLisandro Dalcin PetscFunctionBegin; 461de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 462698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 463*0598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 464*0598e1a2SLisandro Dalcin ierr = GmshNodesCreate(num, &nodes);CHKERRQ(ierr); 465*0598e1a2SLisandro Dalcin mesh->numNodes = num; 466*0598e1a2SLisandro Dalcin mesh->nodelist = nodes; 467*0598e1a2SLisandro Dalcin for (n = 0; n < num; ++n) { 468*0598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 46911c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 47011c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 471*0598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 47211c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 473*0598e1a2SLisandro Dalcin nodes->id[n] = nid; 47411c56cb0SLisandro Dalcin } 47511c56cb0SLisandro Dalcin PetscFunctionReturn(0); 47611c56cb0SLisandro Dalcin } 47711c56cb0SLisandro Dalcin 478de78e4feSLisandro Dalcin /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 479de78e4feSLisandro Dalcin file contents multiple times to figure out the true number of cells and facets 480de78e4feSLisandro Dalcin in the given mesh. To make this more efficient we read the file contents only 481de78e4feSLisandro Dalcin once and store them in memory, while determining the true number of cells. */ 482*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v22(GmshFile* gmsh, GmshMesh *mesh) 48311c56cb0SLisandro Dalcin { 484698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 485698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 486698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 487de78e4feSLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 488*0598e1a2SLisandro Dalcin int i, c, p, num, ibuf[1+4+1000], snum; 489*0598e1a2SLisandro Dalcin int cellType, numElem, numVerts, numNodes, numTags; 490df032b82SMatthew G. Knepley GmshElement *elements; 491*0598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 492df032b82SMatthew G. Knepley PetscErrorCode ierr; 493df032b82SMatthew G. Knepley 494df032b82SMatthew G. Knepley PetscFunctionBegin; 495feb237baSPierre Jolivet ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 496698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &num); 497*0598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 498*0598e1a2SLisandro Dalcin ierr = GmshElementsCreate(num, &elements);CHKERRQ(ierr); 499*0598e1a2SLisandro Dalcin mesh->numElems = num; 500*0598e1a2SLisandro Dalcin mesh->elements = elements; 501698a79b9SLisandro Dalcin for (c = 0; c < num;) { 50211c56cb0SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 50311c56cb0SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 504*0598e1a2SLisandro Dalcin 505*0598e1a2SLisandro Dalcin cellType = binary ? ibuf[0] : ibuf[1]; 506*0598e1a2SLisandro Dalcin numElem = binary ? ibuf[1] : 1; 507df032b82SMatthew G. Knepley numTags = ibuf[2]; 508*0598e1a2SLisandro Dalcin 509*0598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 510*0598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 511*0598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 512*0598e1a2SLisandro Dalcin 513df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 514*0598e1a2SLisandro Dalcin GmshElement *element = elements + c; 515*0598e1a2SLisandro Dalcin const int off = binary ? 0 : 1, nint = 1 + numTags + numNodes - off; 516*0598e1a2SLisandro Dalcin const int *id = ibuf, *nodes = ibuf + 1 + numTags, *tags = ibuf + 1; 517*0598e1a2SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf+off, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 518*0598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf+off, PETSC_ENUM, nint);CHKERRQ(ierr);} 519*0598e1a2SLisandro Dalcin element->id = id[0]; 520*0598e1a2SLisandro Dalcin element->dim = GmshCellMap[cellType].dim; 521*0598e1a2SLisandro Dalcin element->cellType = cellType; 522*0598e1a2SLisandro Dalcin element->numVerts = numVerts; 523*0598e1a2SLisandro Dalcin element->numNodes = numNodes; 524*0598e1a2SLisandro Dalcin element->numTags = PetscMin(numTags, 4); 525*0598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 526*0598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 527*0598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = tags[p]; 528df032b82SMatthew G. Knepley } 529df032b82SMatthew G. Knepley } 530df032b82SMatthew G. Knepley PetscFunctionReturn(0); 531df032b82SMatthew G. Knepley } 5327d282ae0SMichael Lange 533de78e4feSLisandro Dalcin /* 534de78e4feSLisandro Dalcin $Entities 535de78e4feSLisandro Dalcin numPoints(unsigned long) numCurves(unsigned long) 536de78e4feSLisandro Dalcin numSurfaces(unsigned long) numVolumes(unsigned long) 537de78e4feSLisandro Dalcin // points 538de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 539de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 540de78e4feSLisandro Dalcin numPhysicals(unsigned long) phyisicalTag[...](int) 541de78e4feSLisandro Dalcin ... 542de78e4feSLisandro Dalcin // curves 543de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 544de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 545de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 546de78e4feSLisandro Dalcin numBREPVert(unsigned long) tagBREPVert[...](int) 547de78e4feSLisandro Dalcin ... 548de78e4feSLisandro Dalcin // surfaces 549de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 550de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 551de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 552de78e4feSLisandro Dalcin numBREPCurve(unsigned long) tagBREPCurve[...](int) 553de78e4feSLisandro Dalcin ... 554de78e4feSLisandro Dalcin // volumes 555de78e4feSLisandro Dalcin tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 556de78e4feSLisandro Dalcin boxMaxX(double) boxMaxY(double) boxMaxZ(double) 557de78e4feSLisandro Dalcin numPhysicals(unsigned long) physicalTag[...](int) 558de78e4feSLisandro Dalcin numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 559de78e4feSLisandro Dalcin ... 560de78e4feSLisandro Dalcin $EndEntities 561de78e4feSLisandro Dalcin */ 562*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v40(GmshFile *gmsh, GmshMesh *mesh) 563de78e4feSLisandro Dalcin { 564698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 565698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 566698a79b9SLisandro Dalcin long index, num, lbuf[4]; 567730356f6SLisandro Dalcin int dim, eid, numTags, *ibuf, t; 568698a79b9SLisandro Dalcin PetscInt count[4], i; 569a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 570de78e4feSLisandro Dalcin PetscErrorCode ierr; 571de78e4feSLisandro Dalcin 572de78e4feSLisandro Dalcin PetscFunctionBegin; 573698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, lbuf, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 574698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(lbuf, PETSC_LONG, 4);CHKERRQ(ierr);} 575698a79b9SLisandro Dalcin for (i = 0; i < 4; ++i) count[i] = lbuf[i]; 576*0598e1a2SLisandro Dalcin ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr); 577de78e4feSLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 578730356f6SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 579730356f6SLisandro Dalcin ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 580730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 581*0598e1a2SLisandro Dalcin ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 582de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 583de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 584de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 585de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 586698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 587de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 588730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 589de78e4feSLisandro Dalcin entity->numTags = numTags = (int) PetscMin(num, 4); 590de78e4feSLisandro Dalcin for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 591de78e4feSLisandro Dalcin if (dim == 0) continue; 592de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 593de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 594698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, num, sizeof(int), &ibuf);CHKERRQ(ierr); 595de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 596730356f6SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 597de78e4feSLisandro Dalcin } 598de78e4feSLisandro Dalcin } 599de78e4feSLisandro Dalcin PetscFunctionReturn(0); 600de78e4feSLisandro Dalcin } 601de78e4feSLisandro Dalcin 602de78e4feSLisandro Dalcin /* 603de78e4feSLisandro Dalcin $Nodes 604de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numNodes(unsigned long) 605de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 606de78e4feSLisandro Dalcin tag(int) x(double) y(double) z(double) 607de78e4feSLisandro Dalcin ... 608de78e4feSLisandro Dalcin ... 609de78e4feSLisandro Dalcin $EndNodes 610de78e4feSLisandro Dalcin */ 611*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v40(GmshFile *gmsh, GmshMesh *mesh) 612de78e4feSLisandro Dalcin { 613698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 614698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 615*0598e1a2SLisandro Dalcin long block, node, n, numEntityBlocks, numTotalNodes, numNodes; 616de78e4feSLisandro Dalcin int info[3], nid; 617*0598e1a2SLisandro Dalcin GmshNodes *nodes; 618de78e4feSLisandro Dalcin PetscErrorCode ierr; 619de78e4feSLisandro Dalcin 620de78e4feSLisandro Dalcin PetscFunctionBegin; 621de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 622de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 623de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 624de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 625*0598e1a2SLisandro Dalcin ierr = GmshNodesCreate(numTotalNodes, &nodes);CHKERRQ(ierr); 626*0598e1a2SLisandro Dalcin mesh->numNodes = numTotalNodes; 627*0598e1a2SLisandro Dalcin mesh->nodelist = nodes; 628*0598e1a2SLisandro Dalcin for (n = 0, block = 0; block < numEntityBlocks; ++block) { 629de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 630de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 631de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 632698a79b9SLisandro Dalcin if (gmsh->binary) { 633277f51e8SBarry Smith size_t nbytes = sizeof(int) + 3*sizeof(double); 634277f51e8SBarry Smith char *cbuf = NULL; /* dummy value to prevent warning from compiler about possible unitilized value */ 635698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 636de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 637*0598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 638de78e4feSLisandro Dalcin char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 639*0598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 64030815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr);} 64130815ce0SLisandro Dalcin if (!PetscBinaryBigEndian()) {ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 642de78e4feSLisandro Dalcin ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 643de78e4feSLisandro Dalcin ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 644de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 645de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 646*0598e1a2SLisandro Dalcin nodes->id[n] = nid; 647de78e4feSLisandro Dalcin } 648de78e4feSLisandro Dalcin } else { 649*0598e1a2SLisandro Dalcin for (node = 0; node < numNodes; ++node, ++n) { 650*0598e1a2SLisandro Dalcin double *xyz = nodes->xyz + n*3; 651de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 652de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 653*0598e1a2SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 654de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 655*0598e1a2SLisandro Dalcin nodes->id[n] = nid; 656de78e4feSLisandro Dalcin } 657de78e4feSLisandro Dalcin } 658de78e4feSLisandro Dalcin } 659de78e4feSLisandro Dalcin PetscFunctionReturn(0); 660de78e4feSLisandro Dalcin } 661de78e4feSLisandro Dalcin 662de78e4feSLisandro Dalcin /* 663de78e4feSLisandro Dalcin $Elements 664de78e4feSLisandro Dalcin numEntityBlocks(unsigned long) numElements(unsigned long) 665de78e4feSLisandro Dalcin tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 666de78e4feSLisandro Dalcin tag(int) numVert[...](int) 667de78e4feSLisandro Dalcin ... 668de78e4feSLisandro Dalcin ... 669de78e4feSLisandro Dalcin $EndElements 670de78e4feSLisandro Dalcin */ 671*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v40(GmshFile *gmsh, GmshMesh *mesh) 672de78e4feSLisandro Dalcin { 673698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 674698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 675de78e4feSLisandro Dalcin long c, block, numEntityBlocks, numTotalElements, elem, numElements; 676de78e4feSLisandro Dalcin int p, info[3], *ibuf = NULL; 677*0598e1a2SLisandro Dalcin int eid, dim, cellType, numVerts, numNodes, numTags; 678a5ba3d71SLisandro Dalcin GmshEntity *entity = NULL; 679de78e4feSLisandro Dalcin GmshElement *elements; 680*0598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 681de78e4feSLisandro Dalcin PetscErrorCode ierr; 682de78e4feSLisandro Dalcin 683de78e4feSLisandro Dalcin PetscFunctionBegin; 684de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 685de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 686de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 687de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 688*0598e1a2SLisandro Dalcin ierr = GmshElementsCreate(numTotalElements, &elements);CHKERRQ(ierr); 689*0598e1a2SLisandro Dalcin mesh->numElems = numTotalElements; 690*0598e1a2SLisandro Dalcin mesh->elements = elements; 691de78e4feSLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 692de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 693de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 694de78e4feSLisandro Dalcin eid = info[0]; dim = info[1]; cellType = info[2]; 695*0598e1a2SLisandro Dalcin ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr); 696*0598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 697*0598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 698*0598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 699730356f6SLisandro Dalcin numTags = entity->numTags; 700de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 701de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 702698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 703de78e4feSLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 704de78e4feSLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 705de78e4feSLisandro Dalcin for (elem = 0; elem < numElements; ++elem, ++c) { 706de78e4feSLisandro Dalcin GmshElement *element = elements + c; 707*0598e1a2SLisandro Dalcin const int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 708*0598e1a2SLisandro Dalcin element->id = id[0]; 709de78e4feSLisandro Dalcin element->dim = dim; 710de78e4feSLisandro Dalcin element->cellType = cellType; 711*0598e1a2SLisandro Dalcin element->numVerts = numVerts; 712de78e4feSLisandro Dalcin element->numNodes = numNodes; 713de78e4feSLisandro Dalcin element->numTags = numTags; 714*0598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 715*0598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 716*0598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 717de78e4feSLisandro Dalcin } 718de78e4feSLisandro Dalcin } 719698a79b9SLisandro Dalcin PetscFunctionReturn(0); 720698a79b9SLisandro Dalcin } 721698a79b9SLisandro Dalcin 722*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v40(GmshFile *gmsh, PetscInt periodicMap[]) 723698a79b9SLisandro Dalcin { 724698a79b9SLisandro Dalcin PetscViewer viewer = gmsh->viewer; 725698a79b9SLisandro Dalcin int fileFormat = gmsh->fileFormat; 726698a79b9SLisandro Dalcin PetscBool binary = gmsh->binary; 727698a79b9SLisandro Dalcin PetscBool byteSwap = gmsh->byteSwap; 728698a79b9SLisandro Dalcin int numPeriodic, snum, i; 729698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 730*0598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 731698a79b9SLisandro Dalcin PetscErrorCode ierr; 732698a79b9SLisandro Dalcin 733698a79b9SLisandro Dalcin PetscFunctionBegin; 734698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 735698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 736698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numPeriodic); 737*0598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 738698a79b9SLisandro Dalcin } else { 739698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 740698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 741698a79b9SLisandro Dalcin } 742698a79b9SLisandro Dalcin for (i = 0; i < numPeriodic; i++) { 743698a79b9SLisandro Dalcin int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 744698a79b9SLisandro Dalcin long j, nNodes; 745698a79b9SLisandro Dalcin double affine[16]; 746698a79b9SLisandro Dalcin 747698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 748698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 749698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 750*0598e1a2SLisandro Dalcin if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 751698a79b9SLisandro Dalcin } else { 752698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 753698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 754698a79b9SLisandro Dalcin slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 755698a79b9SLisandro Dalcin } 756698a79b9SLisandro Dalcin (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 757698a79b9SLisandro Dalcin 758698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 759698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 760698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 7614c500f23SPierre Jolivet if (snum != 1) { /* discard transformation and try again */ 762698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 763698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 764698a79b9SLisandro Dalcin snum = sscanf(line, "%ld", &nNodes); 765*0598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 766698a79b9SLisandro Dalcin } 767698a79b9SLisandro Dalcin } else { 768698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 769698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 7704c500f23SPierre Jolivet if (nNodes == -1) { /* discard transformation and try again */ 771698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 772698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 773698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 774698a79b9SLisandro Dalcin } 775698a79b9SLisandro Dalcin } 776698a79b9SLisandro Dalcin 777698a79b9SLisandro Dalcin for (j = 0; j < nNodes; j++) { 778698a79b9SLisandro Dalcin if (fileFormat == 22 || !binary) { 779698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 780698a79b9SLisandro Dalcin snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 781*0598e1a2SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 782698a79b9SLisandro Dalcin } else { 783698a79b9SLisandro Dalcin ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 784698a79b9SLisandro Dalcin if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 785698a79b9SLisandro Dalcin slaveNode = ibuf[0]; masterNode = ibuf[1]; 786698a79b9SLisandro Dalcin } 787*0598e1a2SLisandro Dalcin slaveNode = (int) nodeMap[slaveNode]; 788*0598e1a2SLisandro Dalcin masterNode = (int) nodeMap[masterNode]; 789*0598e1a2SLisandro Dalcin periodicMap[slaveNode] = masterNode; 790698a79b9SLisandro Dalcin } 791698a79b9SLisandro Dalcin } 792698a79b9SLisandro Dalcin PetscFunctionReturn(0); 793698a79b9SLisandro Dalcin } 794698a79b9SLisandro Dalcin 795*0598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 796698a79b9SLisandro Dalcin $Entities 797698a79b9SLisandro Dalcin numPoints(size_t) numCurves(size_t) 798698a79b9SLisandro Dalcin numSurfaces(size_t) numVolumes(size_t) 799698a79b9SLisandro Dalcin pointTag(int) X(double) Y(double) Z(double) 800698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 801698a79b9SLisandro Dalcin ... 802698a79b9SLisandro Dalcin curveTag(int) minX(double) minY(double) minZ(double) 803698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 804698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 805698a79b9SLisandro Dalcin numBoundingPoints(size_t) pointTag(int) ... 806698a79b9SLisandro Dalcin ... 807698a79b9SLisandro Dalcin surfaceTag(int) minX(double) minY(double) minZ(double) 808698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 809698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 810698a79b9SLisandro Dalcin numBoundingCurves(size_t) curveTag(int) ... 811698a79b9SLisandro Dalcin ... 812698a79b9SLisandro Dalcin volumeTag(int) minX(double) minY(double) minZ(double) 813698a79b9SLisandro Dalcin maxX(double) maxY(double) maxZ(double) 814698a79b9SLisandro Dalcin numPhysicalTags(size_t) physicalTag(int) ... 815698a79b9SLisandro Dalcin numBoundngSurfaces(size_t) surfaceTag(int) ... 816698a79b9SLisandro Dalcin ... 817698a79b9SLisandro Dalcin $EndEntities 818698a79b9SLisandro Dalcin */ 819*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities_v41(GmshFile *gmsh, GmshMesh *mesh) 820698a79b9SLisandro Dalcin { 821698a79b9SLisandro Dalcin PetscInt count[4], index, numTags, i; 822698a79b9SLisandro Dalcin int dim, eid, *tags = NULL; 823698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 824698a79b9SLisandro Dalcin PetscErrorCode ierr; 825698a79b9SLisandro Dalcin 826698a79b9SLisandro Dalcin PetscFunctionBegin; 827698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, count, 4);CHKERRQ(ierr); 828*0598e1a2SLisandro Dalcin ierr = GmshEntitiesCreate(count, &mesh->entities);CHKERRQ(ierr); 829698a79b9SLisandro Dalcin for (dim = 0; dim < 4; ++dim) { 830698a79b9SLisandro Dalcin for (index = 0; index < count[dim]; ++index) { 831698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &eid, 1);CHKERRQ(ierr); 832*0598e1a2SLisandro Dalcin ierr = GmshEntitiesAdd(mesh->entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 833698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, entity->bbox, (dim == 0) ? 3 : 6);CHKERRQ(ierr); 834698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 835698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 836698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 837698a79b9SLisandro Dalcin entity->numTags = PetscMin(numTags, 4); 838698a79b9SLisandro Dalcin for (i = 0; i < entity->numTags; ++i) entity->tags[i] = tags[i]; 839698a79b9SLisandro Dalcin if (dim == 0) continue; 840698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numTags, 1);CHKERRQ(ierr); 841698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numTags, sizeof(int), &tags);CHKERRQ(ierr); 842698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, tags, numTags);CHKERRQ(ierr); 843698a79b9SLisandro Dalcin } 844698a79b9SLisandro Dalcin } 845698a79b9SLisandro Dalcin PetscFunctionReturn(0); 846698a79b9SLisandro Dalcin } 847698a79b9SLisandro Dalcin 84833a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 849698a79b9SLisandro Dalcin $Nodes 850698a79b9SLisandro Dalcin numEntityBlocks(size_t) numNodes(size_t) 851698a79b9SLisandro Dalcin minNodeTag(size_t) maxNodeTag(size_t) 852698a79b9SLisandro Dalcin entityDim(int) entityTag(int) parametric(int; 0 or 1) numNodesBlock(size_t) 853698a79b9SLisandro Dalcin nodeTag(size_t) 854698a79b9SLisandro Dalcin ... 855698a79b9SLisandro Dalcin x(double) y(double) z(double) 856698a79b9SLisandro Dalcin < u(double; if parametric and entityDim = 1 or entityDim = 2) > 857698a79b9SLisandro Dalcin < v(double; if parametric and entityDim = 2) > 858698a79b9SLisandro Dalcin ... 859698a79b9SLisandro Dalcin ... 860698a79b9SLisandro Dalcin $EndNodes 861698a79b9SLisandro Dalcin */ 862*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes_v41(GmshFile *gmsh, GmshMesh *mesh) 863698a79b9SLisandro Dalcin { 864698a79b9SLisandro Dalcin int info[3]; 865*0598e1a2SLisandro Dalcin PetscInt sizes[4], numEntityBlocks, numNodes, numNodesBlock = 0, block, node; 866*0598e1a2SLisandro Dalcin GmshNodes *nodes; 867698a79b9SLisandro Dalcin PetscErrorCode ierr; 868698a79b9SLisandro Dalcin 869698a79b9SLisandro Dalcin PetscFunctionBegin; 870698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 871*0598e1a2SLisandro Dalcin numEntityBlocks = sizes[0]; numNodes = sizes[1]; 872*0598e1a2SLisandro Dalcin ierr = GmshNodesCreate(numNodes, &nodes);CHKERRQ(ierr); 873*0598e1a2SLisandro Dalcin mesh->numNodes = numNodes; 874*0598e1a2SLisandro Dalcin mesh->nodelist = nodes; 875698a79b9SLisandro Dalcin for (block = 0, node = 0; block < numEntityBlocks; ++block, node += numNodesBlock) { 876698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 877698a79b9SLisandro Dalcin if (info[2] != 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Parametric coordinates not supported"); 878*0598e1a2SLisandro Dalcin ierr = GmshReadSize(gmsh, &numNodesBlock, 1);CHKERRQ(ierr); 879*0598e1a2SLisandro Dalcin ierr = GmshReadSize(gmsh, nodes->id+node, numNodesBlock);CHKERRQ(ierr); 880*0598e1a2SLisandro Dalcin ierr = GmshReadDouble(gmsh, nodes->xyz+node*3, numNodesBlock*3);CHKERRQ(ierr); 881698a79b9SLisandro Dalcin } 882*0598e1a2SLisandro Dalcin gmsh->nodeStart = sizes[2]; 883*0598e1a2SLisandro Dalcin gmsh->nodeEnd = sizes[3]+1; 884698a79b9SLisandro Dalcin PetscFunctionReturn(0); 885698a79b9SLisandro Dalcin } 886698a79b9SLisandro Dalcin 88733a088b6SMatthew G. Knepley /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 888698a79b9SLisandro Dalcin $Elements 889698a79b9SLisandro Dalcin numEntityBlocks(size_t) numElements(size_t) 890698a79b9SLisandro Dalcin minElementTag(size_t) maxElementTag(size_t) 891698a79b9SLisandro Dalcin entityDim(int) entityTag(int) elementType(int; see below) numElementsBlock(size_t) 892698a79b9SLisandro Dalcin elementTag(size_t) nodeTag(size_t) ... 893698a79b9SLisandro Dalcin ... 894698a79b9SLisandro Dalcin ... 895698a79b9SLisandro Dalcin $EndElements 896698a79b9SLisandro Dalcin */ 897*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements_v41(GmshFile *gmsh, GmshMesh *mesh) 898698a79b9SLisandro Dalcin { 899*0598e1a2SLisandro Dalcin int info[3], eid, dim, cellType; 900*0598e1a2SLisandro Dalcin PetscInt sizes[4], *ibuf = NULL, numEntityBlocks, numElements, numBlockElements, numVerts, numNodes, numTags, block, elem, c, p; 901698a79b9SLisandro Dalcin GmshEntity *entity = NULL; 902698a79b9SLisandro Dalcin GmshElement *elements; 903*0598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 904698a79b9SLisandro Dalcin PetscErrorCode ierr; 905698a79b9SLisandro Dalcin 906698a79b9SLisandro Dalcin PetscFunctionBegin; 907698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, sizes, 4);CHKERRQ(ierr); 908698a79b9SLisandro Dalcin numEntityBlocks = sizes[0]; numElements = sizes[1]; 909*0598e1a2SLisandro Dalcin ierr = GmshElementsCreate(numElements, &elements);CHKERRQ(ierr); 910*0598e1a2SLisandro Dalcin mesh->numElems = numElements; 911*0598e1a2SLisandro Dalcin mesh->elements = elements; 912698a79b9SLisandro Dalcin for (c = 0, block = 0; block < numEntityBlocks; ++block) { 913698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 914698a79b9SLisandro Dalcin dim = info[0]; eid = info[1]; cellType = info[2]; 915*0598e1a2SLisandro Dalcin ierr = GmshEntitiesGet(mesh->entities, dim, eid, &entity);CHKERRQ(ierr); 916*0598e1a2SLisandro Dalcin ierr = GmshCellTypeCheck(cellType);CHKERRQ(ierr); 917*0598e1a2SLisandro Dalcin numVerts = GmshCellMap[cellType].numVerts; 918*0598e1a2SLisandro Dalcin numNodes = GmshCellMap[cellType].numNodes; 919698a79b9SLisandro Dalcin numTags = entity->numTags; 920698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numBlockElements, 1);CHKERRQ(ierr); 921698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, (1+numNodes)*numBlockElements, sizeof(PetscInt), &ibuf);CHKERRQ(ierr); 922698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, ibuf, (1+numNodes)*numBlockElements);CHKERRQ(ierr); 923698a79b9SLisandro Dalcin for (elem = 0; elem < numBlockElements; ++elem, ++c) { 924698a79b9SLisandro Dalcin GmshElement *element = elements + c; 925*0598e1a2SLisandro Dalcin const PetscInt *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 926698a79b9SLisandro Dalcin element->id = id[0]; 927698a79b9SLisandro Dalcin element->dim = dim; 928698a79b9SLisandro Dalcin element->cellType = cellType; 929*0598e1a2SLisandro Dalcin element->numVerts = numVerts; 930698a79b9SLisandro Dalcin element->numNodes = numNodes; 931698a79b9SLisandro Dalcin element->numTags = numTags; 932*0598e1a2SLisandro Dalcin ierr = PetscSegBufferGet(mesh->segbuf, (size_t)element->numNodes, &element->nodes);CHKERRQ(ierr); 933*0598e1a2SLisandro Dalcin for (p = 0; p < element->numNodes; p++) element->nodes[p] = nodeMap[nodes[p]]; 934*0598e1a2SLisandro Dalcin for (p = 0; p < element->numTags; p++) element->tags[p] = entity->tags[p]; 935698a79b9SLisandro Dalcin } 936698a79b9SLisandro Dalcin } 937698a79b9SLisandro Dalcin PetscFunctionReturn(0); 938698a79b9SLisandro Dalcin } 939698a79b9SLisandro Dalcin 940*0598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 941698a79b9SLisandro Dalcin $Periodic 942698a79b9SLisandro Dalcin numPeriodicLinks(size_t) 943698a79b9SLisandro Dalcin entityDim(int) entityTag(int) entityTagMaster(int) 944698a79b9SLisandro Dalcin numAffine(size_t) value(double) ... 945698a79b9SLisandro Dalcin numCorrespondingNodes(size_t) 946698a79b9SLisandro Dalcin nodeTag(size_t) nodeTagMaster(size_t) 947698a79b9SLisandro Dalcin ... 948698a79b9SLisandro Dalcin ... 949698a79b9SLisandro Dalcin $EndPeriodic 950698a79b9SLisandro Dalcin */ 951*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic_v41(GmshFile *gmsh, PetscInt periodicMap[]) 952698a79b9SLisandro Dalcin { 953698a79b9SLisandro Dalcin int info[3]; 954698a79b9SLisandro Dalcin double dbuf[16]; 955*0598e1a2SLisandro Dalcin PetscInt numPeriodicLinks, numAffine, numCorrespondingNodes, *nodeTags = NULL, link, node; 956*0598e1a2SLisandro Dalcin PetscInt *nodeMap = gmsh->nodeMap; 957698a79b9SLisandro Dalcin PetscErrorCode ierr; 958698a79b9SLisandro Dalcin 959698a79b9SLisandro Dalcin PetscFunctionBegin; 960698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numPeriodicLinks, 1);CHKERRQ(ierr); 961698a79b9SLisandro Dalcin for (link = 0; link < numPeriodicLinks; ++link) { 962698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, info, 3);CHKERRQ(ierr); 963698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numAffine, 1);CHKERRQ(ierr); 964698a79b9SLisandro Dalcin ierr = GmshReadDouble(gmsh, dbuf, numAffine);CHKERRQ(ierr); 965698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, &numCorrespondingNodes, 1);CHKERRQ(ierr); 966698a79b9SLisandro Dalcin ierr = GmshBufferGet(gmsh, numCorrespondingNodes, sizeof(PetscInt), &nodeTags);CHKERRQ(ierr); 967698a79b9SLisandro Dalcin ierr = GmshReadSize(gmsh, nodeTags, numCorrespondingNodes*2);CHKERRQ(ierr); 968698a79b9SLisandro Dalcin for (node = 0; node < numCorrespondingNodes; ++node) { 969*0598e1a2SLisandro Dalcin PetscInt slaveNode = nodeMap[nodeTags[node*2+0]]; 970*0598e1a2SLisandro Dalcin PetscInt masterNode = nodeMap[nodeTags[node*2+1]]; 971*0598e1a2SLisandro Dalcin periodicMap[slaveNode] = masterNode; 972698a79b9SLisandro Dalcin } 973698a79b9SLisandro Dalcin } 974698a79b9SLisandro Dalcin PetscFunctionReturn(0); 975698a79b9SLisandro Dalcin } 976698a79b9SLisandro Dalcin 977*0598e1a2SLisandro Dalcin /* http://gmsh.info/dev/doc/texinfo/gmsh.html#MSH-file-format 978d6689ca9SLisandro Dalcin $MeshFormat // same as MSH version 2 979d6689ca9SLisandro Dalcin version(ASCII double; currently 4.1) 980d6689ca9SLisandro Dalcin file-type(ASCII int; 0 for ASCII mode, 1 for binary mode) 981d6689ca9SLisandro Dalcin data-size(ASCII int; sizeof(size_t)) 982d6689ca9SLisandro Dalcin < int with value one; only in binary mode, to detect endianness > 983d6689ca9SLisandro Dalcin $EndMeshFormat 984d6689ca9SLisandro Dalcin */ 985*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadMeshFormat(GmshFile *gmsh) 986698a79b9SLisandro Dalcin { 987698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 988698a79b9SLisandro Dalcin int snum, fileType, fileFormat, dataSize, checkEndian; 989698a79b9SLisandro Dalcin float version; 990698a79b9SLisandro Dalcin PetscErrorCode ierr; 991698a79b9SLisandro Dalcin 992698a79b9SLisandro Dalcin PetscFunctionBegin; 993698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 3);CHKERRQ(ierr); 994698a79b9SLisandro Dalcin snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 995*0598e1a2SLisandro Dalcin if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 996*0598e1a2SLisandro 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); 997*0598e1a2SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 998*0598e1a2SLisandro 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); 999*0598e1a2SLisandro Dalcin if (gmsh->binary && !fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is binary but Gmsh file is ASCII"); 1000*0598e1a2SLisandro Dalcin if (!gmsh->binary && fileType) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Viewer is ASCII but Gmsh file is binary"); 1001698a79b9SLisandro Dalcin fileFormat = (int)(version*10); /* XXX Should use (int)roundf(version*10) ? */ 1002*0598e1a2SLisandro 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); 1003*0598e1a2SLisandro 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); 1004698a79b9SLisandro Dalcin gmsh->fileFormat = fileFormat; 1005698a79b9SLisandro Dalcin gmsh->dataSize = dataSize; 1006698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_FALSE; 1007698a79b9SLisandro Dalcin if (gmsh->binary) { 1008698a79b9SLisandro Dalcin ierr = GmshReadInt(gmsh, &checkEndian, 1);CHKERRQ(ierr); 1009698a79b9SLisandro Dalcin if (checkEndian != 1) { 1010698a79b9SLisandro Dalcin ierr = PetscByteSwap(&checkEndian, PETSC_ENUM, 1);CHKERRQ(ierr); 1011*0598e1a2SLisandro Dalcin if (checkEndian != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to detect endianness in Gmsh file header: %s", line); 1012698a79b9SLisandro Dalcin gmsh->byteSwap = PETSC_TRUE; 1013698a79b9SLisandro Dalcin } 1014698a79b9SLisandro Dalcin } 1015698a79b9SLisandro Dalcin PetscFunctionReturn(0); 1016698a79b9SLisandro Dalcin } 1017698a79b9SLisandro Dalcin 10181f643af3SLisandro Dalcin /* 10191f643af3SLisandro Dalcin PhysicalNames 10201f643af3SLisandro Dalcin numPhysicalNames(ASCII int) 10211f643af3SLisandro Dalcin dimension(ASCII int) physicalTag(ASCII int) "name"(127 characters max) 10221f643af3SLisandro Dalcin ... 10231f643af3SLisandro Dalcin $EndPhysicalNames 10241f643af3SLisandro Dalcin */ 1025*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPhysicalNames(GmshFile *gmsh) 1026698a79b9SLisandro Dalcin { 10271f643af3SLisandro Dalcin char line[PETSC_MAX_PATH_LEN], name[128+2], *p, *q; 10281f643af3SLisandro Dalcin int snum, numRegions, region, dim, tag; 1029698a79b9SLisandro Dalcin PetscErrorCode ierr; 1030698a79b9SLisandro Dalcin 1031698a79b9SLisandro Dalcin PetscFunctionBegin; 1032698a79b9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 1);CHKERRQ(ierr); 1033698a79b9SLisandro Dalcin snum = sscanf(line, "%d", &numRegions); 1034*0598e1a2SLisandro Dalcin if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 1035698a79b9SLisandro Dalcin for (region = 0; region < numRegions; ++region) { 10361f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 10371f643af3SLisandro Dalcin snum = sscanf(line, "%d %d", &dim, &tag); 1038*0598e1a2SLisandro Dalcin if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10391f643af3SLisandro Dalcin ierr = GmshReadString(gmsh, line, -(PetscInt)sizeof(line));CHKERRQ(ierr); 10401f643af3SLisandro Dalcin ierr = PetscStrchr(line, '"', &p);CHKERRQ(ierr); 1041*0598e1a2SLisandro Dalcin if (!p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10421f643af3SLisandro Dalcin ierr = PetscStrrchr(line, '"', &q);CHKERRQ(ierr); 1043*0598e1a2SLisandro Dalcin if (q == p) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "File is not a valid Gmsh file"); 10441f643af3SLisandro Dalcin ierr = PetscStrncpy(name, p+1, (size_t)(q-p-1));CHKERRQ(ierr); 1045698a79b9SLisandro Dalcin } 1046de78e4feSLisandro Dalcin PetscFunctionReturn(0); 1047de78e4feSLisandro Dalcin } 1048de78e4feSLisandro Dalcin 1049*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadEntities(GmshFile *gmsh, GmshMesh *mesh) 105096ca5757SLisandro Dalcin { 1051*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 1052*0598e1a2SLisandro Dalcin 1053*0598e1a2SLisandro Dalcin PetscFunctionBegin; 1054*0598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 1055*0598e1a2SLisandro Dalcin case 41: ierr = GmshReadEntities_v41(gmsh, mesh);CHKERRQ(ierr); break; 1056*0598e1a2SLisandro Dalcin default: ierr = GmshReadEntities_v40(gmsh, mesh);CHKERRQ(ierr); break; 105796ca5757SLisandro Dalcin } 1058*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1059*0598e1a2SLisandro Dalcin } 1060*0598e1a2SLisandro Dalcin 1061*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadNodes(GmshFile *gmsh, GmshMesh *mesh) 1062*0598e1a2SLisandro Dalcin { 1063*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 1064*0598e1a2SLisandro Dalcin 1065*0598e1a2SLisandro Dalcin PetscFunctionBegin; 1066*0598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 1067*0598e1a2SLisandro Dalcin case 41: ierr = GmshReadNodes_v41(gmsh, mesh);CHKERRQ(ierr); break; 1068*0598e1a2SLisandro Dalcin case 40: ierr = GmshReadNodes_v40(gmsh, mesh);CHKERRQ(ierr); break; 1069*0598e1a2SLisandro Dalcin default: ierr = GmshReadNodes_v22(gmsh, mesh);CHKERRQ(ierr); break; 1070*0598e1a2SLisandro Dalcin } 1071*0598e1a2SLisandro Dalcin 1072*0598e1a2SLisandro Dalcin { /* Gmsh v2.2/v4.0 does not provide min/max node tags */ 1073*0598e1a2SLisandro Dalcin if (mesh->numNodes > 0 && gmsh->nodeEnd >= gmsh->nodeStart) { 1074*0598e1a2SLisandro Dalcin PetscInt tagMin = PETSC_MAX_INT, tagMax = PETSC_MIN_INT, n; 1075*0598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 1076*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 1077*0598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 1078*0598e1a2SLisandro Dalcin tagMin = PetscMin(tag, tagMin); 1079*0598e1a2SLisandro Dalcin tagMax = PetscMax(tag, tagMax); 1080*0598e1a2SLisandro Dalcin } 1081*0598e1a2SLisandro Dalcin gmsh->nodeStart = tagMin; 1082*0598e1a2SLisandro Dalcin gmsh->nodeEnd = tagMax+1; 1083*0598e1a2SLisandro Dalcin } 1084*0598e1a2SLisandro Dalcin } 1085*0598e1a2SLisandro Dalcin 1086*0598e1a2SLisandro Dalcin { /* Support for sparse node tags */ 1087*0598e1a2SLisandro Dalcin PetscInt n, t; 1088*0598e1a2SLisandro Dalcin GmshNodes *nodes = mesh->nodelist; 1089*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(gmsh->nodeEnd - gmsh->nodeStart, &gmsh->nbuf);CHKERRQ(ierr); 1090*0598e1a2SLisandro Dalcin for (t = 0; t < gmsh->nodeEnd - gmsh->nodeStart; ++t) gmsh->nbuf[t] = PETSC_MIN_INT; 1091*0598e1a2SLisandro Dalcin gmsh->nodeMap = gmsh->nbuf - gmsh->nodeStart; 1092*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) { 1093*0598e1a2SLisandro Dalcin const PetscInt tag = nodes->id[n]; 1094*0598e1a2SLisandro Dalcin if (gmsh->nodeMap[tag] >= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Repeated node tag %D", tag); 1095*0598e1a2SLisandro Dalcin gmsh->nodeMap[tag] = n; 1096*0598e1a2SLisandro Dalcin } 1097*0598e1a2SLisandro Dalcin } 1098*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1099*0598e1a2SLisandro Dalcin } 1100*0598e1a2SLisandro Dalcin 1101*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadElements(GmshFile *gmsh, GmshMesh *mesh) 1102*0598e1a2SLisandro Dalcin { 1103*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 1104*0598e1a2SLisandro Dalcin 1105*0598e1a2SLisandro Dalcin PetscFunctionBegin; 1106*0598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 1107*0598e1a2SLisandro Dalcin case 41: ierr = GmshReadElements_v41(gmsh, mesh);CHKERRQ(ierr); break; 1108*0598e1a2SLisandro Dalcin case 40: ierr = GmshReadElements_v40(gmsh, mesh);CHKERRQ(ierr); break; 1109*0598e1a2SLisandro Dalcin default: ierr = GmshReadElements_v22(gmsh, mesh);CHKERRQ(ierr); break; 1110*0598e1a2SLisandro Dalcin } 1111*0598e1a2SLisandro Dalcin 1112*0598e1a2SLisandro Dalcin { /* Reorder elements by codimension and polytope type */ 1113*0598e1a2SLisandro Dalcin PetscInt ne = mesh->numElems; 1114*0598e1a2SLisandro Dalcin GmshElement *elements = mesh->elements; 1115*0598e1a2SLisandro Dalcin PetscInt keymap[DM_NUM_POLYTOPES], nk = 0; 1116*0598e1a2SLisandro Dalcin PetscInt offset[DM_NUM_POLYTOPES+1], e, k; 1117*0598e1a2SLisandro Dalcin 1118*0598e1a2SLisandro Dalcin for (k = 0; k < DM_NUM_POLYTOPES; ++k) keymap[k] = PETSC_MIN_INT; 1119*0598e1a2SLisandro Dalcin ierr = PetscMemzero(offset,sizeof(offset));CHKERRQ(ierr); 1120*0598e1a2SLisandro Dalcin 1121*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_TETRAHEDRON] = nk++; 1122*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_HEXAHEDRON] = nk++; 1123*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_TRI_PRISM] = nk++; 1124*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_PYRAMID] = nk++; 1125*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_TRIANGLE] = nk++; 1126*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_QUADRILATERAL] = nk++; 1127*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_SEGMENT] = nk++; 1128*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_VERTEX] = nk++; 1129*0598e1a2SLisandro Dalcin keymap[DM_POLYTOPE_UNKNOWN] = nk++; 1130*0598e1a2SLisandro Dalcin 1131*0598e1a2SLisandro Dalcin ierr = GmshElementsCreate(mesh->numElems, &mesh->elements);CHKERRQ(ierr); 1132*0598e1a2SLisandro Dalcin #define key(eid) keymap[GmshCellMap[elements[(eid)].cellType].polytope] 1133*0598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) offset[1+key(e)]++; 1134*0598e1a2SLisandro Dalcin for (k = 1; k < nk; ++k) offset[k] += offset[k-1]; 1135*0598e1a2SLisandro Dalcin for (e = 0; e < ne; ++e) mesh->elements[offset[key(e)]++] = elements[e]; 1136*0598e1a2SLisandro Dalcin #undef key 1137*0598e1a2SLisandro Dalcin ierr = GmshElementsDestroy(&elements);CHKERRQ(ierr); 1138*0598e1a2SLisandro Dalcin 1139*0598e1a2SLisandro Dalcin mesh->dim = mesh->numElems ? mesh->elements[0].dim : 0; 1140*0598e1a2SLisandro Dalcin } 1141*0598e1a2SLisandro Dalcin 1142*0598e1a2SLisandro Dalcin { 1143*0598e1a2SLisandro Dalcin PetscBT vtx; 1144*0598e1a2SLisandro Dalcin PetscInt dim = mesh->dim, e, n, v; 1145*0598e1a2SLisandro Dalcin 1146*0598e1a2SLisandro Dalcin ierr = PetscBTCreate(mesh->numNodes, &vtx);CHKERRQ(ierr); 1147*0598e1a2SLisandro Dalcin 1148*0598e1a2SLisandro Dalcin /* Compute number of cells and set of vertices */ 1149*0598e1a2SLisandro Dalcin mesh->numCells = 0; 1150*0598e1a2SLisandro Dalcin for (e = 0; e < mesh->numElems; ++e) { 1151*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 1152*0598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) mesh->numCells++; 1153*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; v++) { 1154*0598e1a2SLisandro Dalcin ierr = PetscBTSet(vtx, elem->nodes[v]);CHKERRQ(ierr); 1155*0598e1a2SLisandro Dalcin } 1156*0598e1a2SLisandro Dalcin } 1157*0598e1a2SLisandro Dalcin 1158*0598e1a2SLisandro Dalcin /* Compute numbering for vertices */ 1159*0598e1a2SLisandro Dalcin mesh->numVerts = 0; 1160*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(mesh->numNodes, &mesh->vertexMap);CHKERRQ(ierr); 1161*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 1162*0598e1a2SLisandro Dalcin mesh->vertexMap[n] = PetscBTLookup(vtx, n) ? mesh->numVerts++ : PETSC_MIN_INT; 1163*0598e1a2SLisandro Dalcin 1164*0598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&vtx);CHKERRQ(ierr); 1165*0598e1a2SLisandro Dalcin } 1166*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1167*0598e1a2SLisandro Dalcin } 1168*0598e1a2SLisandro Dalcin 1169*0598e1a2SLisandro Dalcin static PetscErrorCode GmshReadPeriodic(GmshFile *gmsh, GmshMesh *mesh) 1170*0598e1a2SLisandro Dalcin { 1171*0598e1a2SLisandro Dalcin PetscInt n; 1172*0598e1a2SLisandro Dalcin PetscErrorCode ierr; 1173*0598e1a2SLisandro Dalcin 1174*0598e1a2SLisandro Dalcin PetscFunctionBegin; 1175*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(mesh->numNodes, &mesh->periodMap);CHKERRQ(ierr); 1176*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) mesh->periodMap[n] = n; 1177*0598e1a2SLisandro Dalcin switch (gmsh->fileFormat) { 1178*0598e1a2SLisandro Dalcin case 41: ierr = GmshReadPeriodic_v41(gmsh, mesh->periodMap);CHKERRQ(ierr); break; 1179*0598e1a2SLisandro Dalcin default: ierr = GmshReadPeriodic_v40(gmsh, mesh->periodMap);CHKERRQ(ierr); break; 1180*0598e1a2SLisandro Dalcin } 1181*0598e1a2SLisandro Dalcin 1182*0598e1a2SLisandro Dalcin /* Find canonical master nodes */ 1183*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 1184*0598e1a2SLisandro Dalcin while (mesh->periodMap[n] != mesh->periodMap[mesh->periodMap[n]]) 1185*0598e1a2SLisandro Dalcin mesh->periodMap[n] = mesh->periodMap[mesh->periodMap[n]]; 1186*0598e1a2SLisandro Dalcin 1187*0598e1a2SLisandro Dalcin /* Renumber vertices (filter out slaves) */ 1188*0598e1a2SLisandro Dalcin mesh->numVerts = 0; 1189*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 1190*0598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 1191*0598e1a2SLisandro Dalcin if (mesh->periodMap[n] == n) /* is master */ 1192*0598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->numVerts++; 1193*0598e1a2SLisandro Dalcin for (n = 0; n < mesh->numNodes; ++n) 1194*0598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) /* is vertex */ 1195*0598e1a2SLisandro Dalcin if (mesh->periodMap[n] != n) /* is slave */ 1196*0598e1a2SLisandro Dalcin mesh->vertexMap[n] = mesh->vertexMap[mesh->periodMap[n]]; 1197*0598e1a2SLisandro Dalcin PetscFunctionReturn(0); 1198*0598e1a2SLisandro Dalcin } 1199*0598e1a2SLisandro Dalcin 1200*0598e1a2SLisandro Dalcin 1201*0598e1a2SLisandro Dalcin PETSC_STATIC_INLINE DMPolytopeType DMPolytopeTypeFromGmsh(PetscInt cellType) 1202*0598e1a2SLisandro Dalcin { 1203*0598e1a2SLisandro Dalcin return GmshCellMap[cellType].polytope; 120496ca5757SLisandro Dalcin } 120596ca5757SLisandro Dalcin 1206d6689ca9SLisandro Dalcin /*@C 1207d6689ca9SLisandro Dalcin DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 1208d6689ca9SLisandro Dalcin 1209d6689ca9SLisandro Dalcin + comm - The MPI communicator 1210d6689ca9SLisandro Dalcin . filename - Name of the Gmsh file 1211d6689ca9SLisandro Dalcin - interpolate - Create faces and edges in the mesh 1212d6689ca9SLisandro Dalcin 1213d6689ca9SLisandro Dalcin Output Parameter: 1214d6689ca9SLisandro Dalcin . dm - The DM object representing the mesh 1215d6689ca9SLisandro Dalcin 1216d6689ca9SLisandro Dalcin Level: beginner 1217d6689ca9SLisandro Dalcin 1218d6689ca9SLisandro Dalcin .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 1219d6689ca9SLisandro Dalcin @*/ 1220d6689ca9SLisandro Dalcin PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1221d6689ca9SLisandro Dalcin { 1222d6689ca9SLisandro Dalcin PetscViewer viewer; 1223d6689ca9SLisandro Dalcin PetscMPIInt rank; 1224d6689ca9SLisandro Dalcin int fileType; 1225d6689ca9SLisandro Dalcin PetscViewerType vtype; 1226d6689ca9SLisandro Dalcin PetscErrorCode ierr; 1227d6689ca9SLisandro Dalcin 1228d6689ca9SLisandro Dalcin PetscFunctionBegin; 1229d6689ca9SLisandro Dalcin ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1230d6689ca9SLisandro Dalcin 1231d6689ca9SLisandro Dalcin /* Determine Gmsh file type (ASCII or binary) from file header */ 1232d6689ca9SLisandro Dalcin if (!rank) { 1233*0598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1234d6689ca9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1235d6689ca9SLisandro Dalcin int snum; 1236d6689ca9SLisandro Dalcin float version; 1237d6689ca9SLisandro Dalcin 1238580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 1239d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(PETSC_COMM_SELF, &gmsh->viewer);CHKERRQ(ierr); 1240d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(gmsh->viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 1241d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(gmsh->viewer, FILE_MODE_READ);CHKERRQ(ierr); 1242d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(gmsh->viewer, filename);CHKERRQ(ierr); 1243d6689ca9SLisandro Dalcin /* Read only the first two lines of the Gmsh file */ 1244d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1245d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 1246d6689ca9SLisandro Dalcin ierr = GmshReadString(gmsh, line, 2);CHKERRQ(ierr); 1247d6689ca9SLisandro Dalcin snum = sscanf(line, "%f %d", &version, &fileType); 1248*0598e1a2SLisandro Dalcin if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Unable to parse Gmsh file header: %s", line); 1249*0598e1a2SLisandro 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); 1250*0598e1a2SLisandro Dalcin if ((int)version == 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Gmsh file version %3.1f not supported", (double)version); 1251*0598e1a2SLisandro 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); 1252d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&gmsh->viewer);CHKERRQ(ierr); 1253d6689ca9SLisandro Dalcin } 1254d6689ca9SLisandro Dalcin ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 1255d6689ca9SLisandro Dalcin vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 1256d6689ca9SLisandro Dalcin 1257d6689ca9SLisandro Dalcin /* Create appropriate viewer and build plex */ 1258d6689ca9SLisandro Dalcin ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 1259d6689ca9SLisandro Dalcin ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 1260d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 1261d6689ca9SLisandro Dalcin ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 1262d6689ca9SLisandro Dalcin ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 1263d6689ca9SLisandro Dalcin ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1264d6689ca9SLisandro Dalcin PetscFunctionReturn(0); 1265d6689ca9SLisandro Dalcin } 1266d6689ca9SLisandro Dalcin 1267331830f3SMatthew G. Knepley /*@ 12687d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 1269331830f3SMatthew G. Knepley 1270d083f849SBarry Smith Collective 1271331830f3SMatthew G. Knepley 1272331830f3SMatthew G. Knepley Input Parameters: 1273331830f3SMatthew G. Knepley + comm - The MPI communicator 1274331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 1275331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 1276331830f3SMatthew G. Knepley 1277331830f3SMatthew G. Knepley Output Parameter: 1278331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 1279331830f3SMatthew G. Knepley 1280698a79b9SLisandro Dalcin Note: http://gmsh.info/doc/texinfo/gmsh.html#MSH-file-format 1281331830f3SMatthew G. Knepley 1282331830f3SMatthew G. Knepley Level: beginner 1283331830f3SMatthew G. Knepley 1284331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 1285331830f3SMatthew G. Knepley @*/ 1286331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1287331830f3SMatthew G. Knepley { 1288*0598e1a2SLisandro Dalcin GmshMesh *mesh = NULL; 128911c56cb0SLisandro Dalcin PetscViewer parentviewer = NULL; 1290*0598e1a2SLisandro Dalcin PetscBT periodicVerts = NULL; 1291*0598e1a2SLisandro Dalcin PetscBT periodicCells = NULL; 1292331830f3SMatthew G. Knepley PetscSection coordSection; 1293331830f3SMatthew G. Knepley Vec coordinates; 1294*0598e1a2SLisandro Dalcin double *coordsIn; 129584572febSToby Isaac PetscScalar *coords; 1296*0598e1a2SLisandro Dalcin PetscInt dim = 0, coordDim = -1; 1297*0598e1a2SLisandro Dalcin PetscInt numNodes = 0, numElems = 0, numVerts = 0, numCells = 0; 1298*0598e1a2SLisandro Dalcin PetscInt coordSize, *vertexMapInv, cell, cone[8], e, n, v, d; 1299*0598e1a2SLisandro Dalcin PetscBool binary, usemarker = PETSC_FALSE; 1300*0598e1a2SLisandro Dalcin PetscBool hybrid = interpolate, periodic = PETSC_TRUE; 130196ca5757SLisandro Dalcin PetscBool hasTetra = PETSC_FALSE; 1302*0598e1a2SLisandro Dalcin PetscMPIInt rank; 1303331830f3SMatthew G. Knepley PetscErrorCode ierr; 1304331830f3SMatthew G. Knepley 1305331830f3SMatthew G. Knepley PetscFunctionBegin; 1306331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1307ef12879bSLisandro Dalcin ierr = PetscObjectOptionsBegin((PetscObject)viewer);CHKERRQ(ierr); 1308ef12879bSLisandro Dalcin ierr = PetscOptionsHead(PetscOptionsObject,"DMPlex Gmsh options");CHKERRQ(ierr); 1309*0598e1a2SLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_hybrid", "Generate hybrid cell bounds", "DMPlexCreateGmsh", hybrid, &hybrid, NULL);CHKERRQ(ierr); 1310ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_periodic","Read Gmsh periodic section", "DMPlexCreateGmsh", periodic, &periodic, NULL);CHKERRQ(ierr); 1311ef12879bSLisandro Dalcin ierr = PetscOptionsBool("-dm_plex_gmsh_use_marker", "Generate marker label", "DMPlexCreateGmsh", usemarker, &usemarker, NULL);CHKERRQ(ierr); 1312*0598e1a2SLisandro Dalcin ierr = PetscOptionsBoundedInt("-dm_plex_gmsh_spacedim", "Embedding space dimension", "DMPlexCreateGmsh", coordDim, &coordDim, NULL,PETSC_DECIDE);CHKERRQ(ierr); 1313ef12879bSLisandro Dalcin ierr = PetscOptionsTail();CHKERRQ(ierr); 1314ef12879bSLisandro Dalcin ierr = PetscOptionsEnd();CHKERRQ(ierr); 1315*0598e1a2SLisandro Dalcin 1316*0598e1a2SLisandro Dalcin ierr = GmshCellInfoSetUp();CHKERRQ(ierr); 131711c56cb0SLisandro Dalcin 1318331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1319331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1320*0598e1a2SLisandro Dalcin ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr); 132111c56cb0SLisandro Dalcin 132211c56cb0SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 132311c56cb0SLisandro Dalcin 132411c56cb0SLisandro Dalcin /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 13253b46f529SLisandro Dalcin if (binary) { 132611c56cb0SLisandro Dalcin parentviewer = viewer; 132711c56cb0SLisandro Dalcin ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 132811c56cb0SLisandro Dalcin } 132911c56cb0SLisandro Dalcin 13303b46f529SLisandro Dalcin if (!rank) { 1331*0598e1a2SLisandro Dalcin GmshFile gmsh[1]; 1332698a79b9SLisandro Dalcin char line[PETSC_MAX_PATH_LEN]; 1333698a79b9SLisandro Dalcin PetscBool match; 1334331830f3SMatthew G. Knepley 1335580bdb30SBarry Smith ierr = PetscArrayzero(gmsh,1);CHKERRQ(ierr); 1336698a79b9SLisandro Dalcin gmsh->viewer = viewer; 1337698a79b9SLisandro Dalcin gmsh->binary = binary; 1338698a79b9SLisandro Dalcin 1339*0598e1a2SLisandro Dalcin ierr = GmshMeshCreate(&mesh);CHKERRQ(ierr); 1340*0598e1a2SLisandro Dalcin 1341698a79b9SLisandro Dalcin /* Read mesh format */ 1342d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1343d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$MeshFormat", line);CHKERRQ(ierr); 1344*0598e1a2SLisandro Dalcin ierr = GmshReadMeshFormat(gmsh);CHKERRQ(ierr); 1345d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndMeshFormat", line);CHKERRQ(ierr); 134611c56cb0SLisandro Dalcin 1347dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 1348d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1349d6689ca9SLisandro Dalcin ierr = GmshMatch(gmsh, "$PhysicalNames", line, &match);CHKERRQ(ierr); 1350dc0ede3bSMatthew G. Knepley if (match) { 1351*0598e1a2SLisandro Dalcin ierr = GmshExpect(gmsh, "$PhysicalNames", line);CHKERRQ(ierr); 1352*0598e1a2SLisandro Dalcin ierr = GmshReadPhysicalNames(gmsh);CHKERRQ(ierr); 1353d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPhysicalNames", line);CHKERRQ(ierr); 1354698a79b9SLisandro Dalcin /* Initial read for entity section */ 1355d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1356dc0ede3bSMatthew G. Knepley } 135711c56cb0SLisandro Dalcin 1358de78e4feSLisandro Dalcin /* Read entities */ 1359698a79b9SLisandro Dalcin if (gmsh->fileFormat >= 40) { 1360d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Entities", line);CHKERRQ(ierr); 1361*0598e1a2SLisandro Dalcin ierr = GmshReadEntities(gmsh, mesh);CHKERRQ(ierr); 1362d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndEntities", line);CHKERRQ(ierr); 1363698a79b9SLisandro Dalcin /* Initial read for nodes section */ 1364d6689ca9SLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1365de78e4feSLisandro Dalcin } 1366de78e4feSLisandro Dalcin 1367698a79b9SLisandro Dalcin /* Read nodes */ 1368d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Nodes", line);CHKERRQ(ierr); 1369*0598e1a2SLisandro Dalcin ierr = GmshReadNodes(gmsh, mesh);CHKERRQ(ierr); 1370d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndNodes", line);CHKERRQ(ierr); 137111c56cb0SLisandro Dalcin 1372698a79b9SLisandro Dalcin /* Read elements */ 1373feb237baSPierre Jolivet ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1374d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Elements", line);CHKERRQ(ierr); 1375*0598e1a2SLisandro Dalcin ierr = GmshReadElements(gmsh, mesh);CHKERRQ(ierr); 1376d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndElements", line);CHKERRQ(ierr); 1377de78e4feSLisandro Dalcin 1378*0598e1a2SLisandro Dalcin /* Read periodic section (OPTIONAL) */ 1379698a79b9SLisandro Dalcin if (periodic) { 1380ef12879bSLisandro Dalcin ierr = GmshReadSection(gmsh, line);CHKERRQ(ierr); 1381ef12879bSLisandro Dalcin ierr = GmshMatch(gmsh, "$Periodic", line, &periodic);CHKERRQ(ierr); 1382ef12879bSLisandro Dalcin } 1383ef12879bSLisandro Dalcin if (periodic) { 1384d6689ca9SLisandro Dalcin ierr = GmshExpect(gmsh, "$Periodic", line);CHKERRQ(ierr); 1385*0598e1a2SLisandro Dalcin ierr = GmshReadPeriodic(gmsh, mesh);CHKERRQ(ierr); 1386d6689ca9SLisandro Dalcin ierr = GmshReadEndSection(gmsh, "$EndPeriodic", line);CHKERRQ(ierr); 1387698a79b9SLisandro Dalcin } 1388698a79b9SLisandro Dalcin 1389698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->wbuf);CHKERRQ(ierr); 1390698a79b9SLisandro Dalcin ierr = PetscFree(gmsh->sbuf);CHKERRQ(ierr); 1391*0598e1a2SLisandro Dalcin ierr = PetscFree(gmsh->nbuf);CHKERRQ(ierr); 1392*0598e1a2SLisandro Dalcin 1393*0598e1a2SLisandro Dalcin dim = mesh->dim; 1394*0598e1a2SLisandro Dalcin numNodes = mesh->numNodes; 1395*0598e1a2SLisandro Dalcin numElems = mesh->numElems; 1396*0598e1a2SLisandro Dalcin numVerts = mesh->numVerts; 1397*0598e1a2SLisandro Dalcin numCells = mesh->numCells; 1398698a79b9SLisandro Dalcin } 1399698a79b9SLisandro Dalcin 1400698a79b9SLisandro Dalcin if (parentviewer) { 1401698a79b9SLisandro Dalcin ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 1402698a79b9SLisandro Dalcin } 1403698a79b9SLisandro Dalcin 1404*0598e1a2SLisandro Dalcin ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1405*0598e1a2SLisandro Dalcin ierr = MPI_Bcast(&periodic, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 1406698a79b9SLisandro Dalcin 1407*0598e1a2SLisandro Dalcin /* Flag presence of tetrahedra to special case wedges */ 1408*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1409*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1410*0598e1a2SLisandro Dalcin DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType); 1411*0598e1a2SLisandro Dalcin if (ctype == DM_POLYTOPE_TETRAHEDRON) hasTetra = PETSC_TRUE; 1412dea2a3c7SStefano Zampini } 1413dea2a3c7SStefano Zampini 1414*0598e1a2SLisandro Dalcin /* We do not want this label automatically computed, instead we fill it here */ 1415*0598e1a2SLisandro Dalcin ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr); 141611c56cb0SLisandro Dalcin 1417a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 1418*0598e1a2SLisandro Dalcin ierr = DMPlexSetChart(*dm, 0, numCells+numVerts);CHKERRQ(ierr); 1419*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1420*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1421*0598e1a2SLisandro Dalcin DMPolytopeType ctype = DMPolytopeTypeFromGmsh(elem->cellType); 1422*0598e1a2SLisandro Dalcin if (hybrid && hasTetra && ctype == DM_POLYTOPE_TRI_PRISM) ctype = DM_POLYTOPE_TRI_PRISM_TENSOR; 1423*0598e1a2SLisandro Dalcin ierr = DMPlexSetConeSize(*dm, cell, elem->numVerts);CHKERRQ(ierr); 1424*0598e1a2SLisandro Dalcin ierr = DMPlexSetCellType(*dm, cell, ctype);CHKERRQ(ierr); 1425db397164SMichael Lange } 1426*0598e1a2SLisandro Dalcin for (v = numCells; v < numCells+numVerts; ++v) { 142796ca5757SLisandro Dalcin ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr); 142896ca5757SLisandro Dalcin } 1429331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 143096ca5757SLisandro Dalcin 1431a4bb7517SMichael Lange /* Add cell-vertex connections */ 1432*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1433*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1434*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 1435*0598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 1436*0598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1437*0598e1a2SLisandro Dalcin cone[v] = numCells + vv; 1438db397164SMichael Lange } 1439*0598e1a2SLisandro Dalcin ierr = DMPlexReorderCell(*dm, cell, cone);CHKERRQ(ierr); 1440*0598e1a2SLisandro Dalcin ierr = DMPlexSetCone(*dm, cell, cone);CHKERRQ(ierr); 1441a4bb7517SMichael Lange } 144296ca5757SLisandro Dalcin 1443c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1444331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1445331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1446331830f3SMatthew G. Knepley if (interpolate) { 14475fd9971aSMatthew G. Knepley DM idm; 1448331830f3SMatthew G. Knepley 1449331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1450331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 1451331830f3SMatthew G. Knepley *dm = idm; 1452331830f3SMatthew G. Knepley } 1453917266a4SLisandro Dalcin 1454f6c8a31fSLisandro Dalcin if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 1455917266a4SLisandro Dalcin if (!rank && usemarker) { 1456d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 1457d3f73514SStefano Zampini 1458d3f73514SStefano Zampini ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 1459d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 1460d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 1461d3f73514SStefano Zampini PetscInt suppSize; 1462d3f73514SStefano Zampini 1463d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 1464d3f73514SStefano Zampini if (suppSize == 1) { 1465d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 1466d3f73514SStefano Zampini 1467d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1468d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 1469d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 1470d3f73514SStefano Zampini } 1471d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 1472d3f73514SStefano Zampini } 1473d3f73514SStefano Zampini } 1474d3f73514SStefano Zampini } 147516ad7e67SMichael Lange 147616ad7e67SMichael Lange if (!rank) { 147711c56cb0SLisandro Dalcin PetscInt vStart, vEnd; 1478d1a54cefSMatthew G. Knepley 147916ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 1480*0598e1a2SLisandro Dalcin for (cell = 0, e = 0; e < numElems; ++e) { 1481*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + e; 14826e1dd89aSLawrence Mitchell 14836e1dd89aSLawrence Mitchell /* Create cell sets */ 1484*0598e1a2SLisandro Dalcin if (elem->dim == dim && dim > 0) { 1485*0598e1a2SLisandro Dalcin if (elem->numTags > 0) { 1486*0598e1a2SLisandro Dalcin ierr = DMSetLabelValue(*dm, "Cell Sets", cell, elem->tags[0]);CHKERRQ(ierr); 14876e1dd89aSLawrence Mitchell } 1488917266a4SLisandro Dalcin cell++; 14896e1dd89aSLawrence Mitchell } 14900c070f12SSander Arens 1491*0598e1a2SLisandro Dalcin /* Create face sets */ 1492*0598e1a2SLisandro Dalcin if (interpolate && elem->dim == dim-1) { 1493*0598e1a2SLisandro Dalcin PetscInt joinSize; 1494*0598e1a2SLisandro Dalcin const PetscInt *join = NULL; 1495*0598e1a2SLisandro Dalcin /* Find the relevant facet with vertex joins */ 1496*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 1497*0598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[v]; 1498*0598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1499*0598e1a2SLisandro Dalcin cone[v] = vStart + vv; 1500*0598e1a2SLisandro Dalcin } 1501*0598e1a2SLisandro Dalcin ierr = DMPlexGetFullJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr); 1502*0598e1a2SLisandro 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); 1503*0598e1a2SLisandro Dalcin ierr = DMSetLabelValue(*dm, "Face Sets", join[0], elem->tags[0]);CHKERRQ(ierr); 1504*0598e1a2SLisandro Dalcin ierr = DMPlexRestoreJoin(*dm, elem->numVerts, cone, &joinSize, &join);CHKERRQ(ierr); 1505*0598e1a2SLisandro Dalcin } 1506*0598e1a2SLisandro Dalcin 15070c070f12SSander Arens /* Create vertex sets */ 1508*0598e1a2SLisandro Dalcin if (elem->dim == 0) { 1509*0598e1a2SLisandro Dalcin if (elem->numTags > 0) { 1510*0598e1a2SLisandro Dalcin const PetscInt nn = elem->nodes[0]; 1511*0598e1a2SLisandro Dalcin const PetscInt vv = mesh->vertexMap[nn]; 1512*0598e1a2SLisandro Dalcin ierr = DMSetLabelValue(*dm, "Vertex Sets", vStart + vv, elem->tags[0]);CHKERRQ(ierr); 1513*0598e1a2SLisandro Dalcin } 1514*0598e1a2SLisandro Dalcin } 1515*0598e1a2SLisandro Dalcin } 1516*0598e1a2SLisandro Dalcin } 1517*0598e1a2SLisandro Dalcin 1518*0598e1a2SLisandro Dalcin if (periodic) { 1519*0598e1a2SLisandro Dalcin ierr = PetscBTCreate(numVerts, &periodicVerts);CHKERRQ(ierr); 1520*0598e1a2SLisandro Dalcin for (n = 0; n < numNodes; ++n) { 1521*0598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) { 1522*0598e1a2SLisandro Dalcin if (PetscUnlikely(mesh->periodMap[n] != n)) { 1523*0598e1a2SLisandro Dalcin PetscInt m = mesh->periodMap[n]; 1524*0598e1a2SLisandro Dalcin ierr= PetscBTSet(periodicVerts, mesh->vertexMap[n]);CHKERRQ(ierr); 1525*0598e1a2SLisandro Dalcin ierr= PetscBTSet(periodicVerts, mesh->vertexMap[m]);CHKERRQ(ierr); 1526*0598e1a2SLisandro Dalcin } 1527*0598e1a2SLisandro Dalcin } 1528*0598e1a2SLisandro Dalcin } 1529*0598e1a2SLisandro Dalcin ierr = PetscBTCreate(numCells, &periodicCells);CHKERRQ(ierr); 1530*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1531*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1532*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) { 1533*0598e1a2SLisandro Dalcin PetscInt nn = elem->nodes[v]; 1534*0598e1a2SLisandro Dalcin PetscInt vv = mesh->vertexMap[nn]; 1535*0598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicVerts, vv))) { 1536*0598e1a2SLisandro Dalcin ierr = PetscBTSet(periodicCells, cell);CHKERRQ(ierr); break; 15370c070f12SSander Arens } 15380c070f12SSander Arens } 15390c070f12SSander Arens } 154016ad7e67SMichael Lange } 154116ad7e67SMichael Lange 154211c56cb0SLisandro Dalcin /* Create coordinates */ 1543*0598e1a2SLisandro Dalcin if (coordDim < 0) coordDim = dim; 1544*0598e1a2SLisandro Dalcin ierr = DMSetCoordinateDim(*dm, coordDim);CHKERRQ(ierr); 1545979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1546331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1547*0598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldComponents(coordSection, 0, coordDim);CHKERRQ(ierr); 1548*0598e1a2SLisandro Dalcin if (periodic) { /* we need to localize coordinates on cells */ 1549*0598e1a2SLisandro Dalcin ierr = PetscSectionSetChart(coordSection, 0, numCells+numVerts);CHKERRQ(ierr); 1550f45c9589SStefano Zampini } else { 1551*0598e1a2SLisandro Dalcin ierr = PetscSectionSetChart(coordSection, numCells, numCells+numVerts);CHKERRQ(ierr); 1552f45c9589SStefano Zampini } 1553*0598e1a2SLisandro Dalcin if (periodic) { 1554*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1555*0598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 1556*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1557*0598e1a2SLisandro Dalcin PetscInt dof = elem->numVerts * coordDim; 1558*0598e1a2SLisandro Dalcin ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr); 1559*0598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr); 1560f45c9589SStefano Zampini } 1561f45c9589SStefano Zampini } 1562f45c9589SStefano Zampini } 1563*0598e1a2SLisandro Dalcin for (v = numCells; v < numCells+numVerts; ++v) { 1564*0598e1a2SLisandro Dalcin ierr = PetscSectionSetDof(coordSection, v, coordDim);CHKERRQ(ierr); 1565*0598e1a2SLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, v, 0, coordDim);CHKERRQ(ierr); 1566*0598e1a2SLisandro Dalcin } 1567331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1568331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 15698b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1570331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1571331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1572*0598e1a2SLisandro Dalcin ierr = VecSetBlockSize(coordinates, coordDim);CHKERRQ(ierr); 1573331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 1574331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1575*0598e1a2SLisandro Dalcin coordsIn = mesh ? mesh->nodelist->xyz : NULL; 1576*0598e1a2SLisandro Dalcin if (periodic) { 1577*0598e1a2SLisandro Dalcin for (cell = 0; cell < numCells; ++cell) { 1578*0598e1a2SLisandro Dalcin if (PetscUnlikely(PetscBTLookup(periodicCells, cell))) { 1579*0598e1a2SLisandro Dalcin GmshElement *elem = mesh->elements + cell; 1580*0598e1a2SLisandro Dalcin PetscInt off, node; 1581*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 1582*0598e1a2SLisandro Dalcin cone[v] = elem->nodes[v]; 1583*0598e1a2SLisandro Dalcin ierr = DMPlexReorderCell(*dm, cell, cone);CHKERRQ(ierr); 1584*0598e1a2SLisandro Dalcin ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr); 1585*0598e1a2SLisandro Dalcin for (v = 0; v < elem->numVerts; ++v) 1586*0598e1a2SLisandro Dalcin for (node = cone[v], d = 0; d < coordDim; ++d) 1587*0598e1a2SLisandro Dalcin coords[off++] = (PetscReal) coordsIn[node*3+d]; 1588331830f3SMatthew G. Knepley } 1589331830f3SMatthew G. Knepley } 1590331830f3SMatthew G. Knepley } 1591*0598e1a2SLisandro Dalcin ierr = PetscMalloc1(numVerts, &vertexMapInv);CHKERRQ(ierr); 1592*0598e1a2SLisandro Dalcin for (n = 0; n < numNodes; n++) 1593*0598e1a2SLisandro Dalcin if (mesh->vertexMap[n] >= 0) 1594*0598e1a2SLisandro Dalcin vertexMapInv[mesh->vertexMap[n]] = n; 1595*0598e1a2SLisandro Dalcin for (v = 0; v < numVerts; ++v) { 1596*0598e1a2SLisandro Dalcin PetscInt off, node = vertexMapInv[v]; 1597*0598e1a2SLisandro Dalcin ierr = PetscSectionGetOffset(coordSection, numCells + v, &off);CHKERRQ(ierr); 1598*0598e1a2SLisandro Dalcin for (d = 0; d < coordDim; ++d) 1599*0598e1a2SLisandro Dalcin coords[off+d] = (PetscReal) coordsIn[node*3+d]; 1600*0598e1a2SLisandro Dalcin } 1601*0598e1a2SLisandro Dalcin ierr = PetscFree(vertexMapInv);CHKERRQ(ierr); 1602331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1603331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1604*0598e1a2SLisandro Dalcin ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 160511c56cb0SLisandro Dalcin ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 160611c56cb0SLisandro Dalcin 1607*0598e1a2SLisandro Dalcin ierr = GmshMeshDestroy(&mesh);CHKERRQ(ierr); 1608*0598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&periodicVerts);CHKERRQ(ierr); 1609*0598e1a2SLisandro Dalcin ierr = PetscBTDestroy(&periodicCells);CHKERRQ(ierr); 161011c56cb0SLisandro Dalcin 1611*0598e1a2SLisandro Dalcin ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,NULL,NULL,NULL);CHKERRQ(ierr); 1612331830f3SMatthew G. Knepley PetscFunctionReturn(0); 1613331830f3SMatthew G. Knepley } 1614