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