1 #define PETSCDM_DLL 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 #include <petsc/private/hashmapi.h> 4 5 /*@C 6 DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 7 8 + comm - The MPI communicator 9 . filename - Name of the Gmsh file 10 - interpolate - Create faces and edges in the mesh 11 12 Output Parameter: 13 . dm - The DM object representing the mesh 14 15 Level: beginner 16 17 .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 18 @*/ 19 PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 20 { 21 PetscViewer viewer; 22 PetscMPIInt rank; 23 int fileType; 24 PetscViewerType vtype; 25 PetscErrorCode ierr; 26 27 PetscFunctionBegin; 28 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 29 30 /* Determine Gmsh file type (ASCII or binary) from file header */ 31 if (!rank) { 32 PetscViewer vheader; 33 char line[PETSC_MAX_PATH_LEN]; 34 PetscBool match; 35 int snum; 36 float version; 37 38 ierr = PetscViewerCreate(PETSC_COMM_SELF, &vheader);CHKERRQ(ierr); 39 ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr); 40 ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr); 41 ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr); 42 /* Read only the first two lines of the Gmsh file */ 43 ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 44 ierr = PetscStrncmp(line, "$MeshFormat", sizeof(line), &match);CHKERRQ(ierr); 45 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 46 ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 47 snum = sscanf(line, "%f %d", &version, &fileType); 48 if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 49 if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 50 if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported"); 51 if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0"); 52 ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 53 } 54 ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 55 vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 56 57 /* Create appropriate viewer and build plex */ 58 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 59 ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 60 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 61 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 62 ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 63 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 64 PetscFunctionReturn(0); 65 } 66 67 typedef struct { 68 size_t size; 69 void *buffer; 70 } GmshWorkBuffer; 71 72 static PetscErrorCode GmshWorkBufferInit(GmshWorkBuffer *ctx) 73 { 74 PetscFunctionBegin; 75 ctx->size = 0; 76 ctx->buffer = NULL; 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode GmshWorkBufferGet(GmshWorkBuffer *ctx, size_t count, size_t eltsize, void *buf) 81 { 82 size_t size = count*eltsize; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 if (ctx->size < size) { 87 ierr = PetscFree(ctx->buffer);CHKERRQ(ierr); 88 ierr = PetscMalloc(size, &ctx->buffer);CHKERRQ(ierr); 89 ctx->size = size; 90 } 91 *(void**)buf = size ? ctx->buffer : NULL; 92 PetscFunctionReturn(0); 93 } 94 95 static PetscErrorCode GmshWorkBufferFree(GmshWorkBuffer *ctx) 96 { 97 PetscErrorCode ierr; 98 PetscFunctionBegin; 99 ierr = PetscFree(ctx->buffer);CHKERRQ(ierr); 100 PetscFunctionReturn(0); 101 } 102 103 typedef struct { 104 PetscInt dim; /* Entity dimension */ 105 PetscInt id; /* Entity number */ 106 double bbox[6]; /* Bounding box */ 107 PetscInt numTags; /* Size of tag array */ 108 int tags[4]; /* Tag array */ 109 } GmshEntity; 110 111 typedef struct { 112 GmshEntity *entity[4]; 113 PetscHMapI entityMap[4]; 114 } GmshEntities; 115 116 static PetscErrorCode GmshEntitiesCreate(long count[4], GmshEntities **entities) 117 { 118 PetscInt dim; 119 PetscErrorCode ierr; 120 121 PetscFunctionBegin; 122 ierr = PetscNew(entities);CHKERRQ(ierr); 123 for (dim = 0; dim < 4; ++dim) { 124 ierr = PetscCalloc1(count[dim], &(*entities)->entity[dim]);CHKERRQ(ierr); 125 ierr = PetscHMapICreate(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 126 } 127 PetscFunctionReturn(0); 128 } 129 130 static PetscErrorCode GmshEntitiesAdd(GmshEntities *entities, PetscInt index, PetscInt dim, PetscInt eid, GmshEntity** entity) 131 { 132 PetscErrorCode ierr; 133 PetscFunctionBegin; 134 ierr = PetscHMapISet(entities->entityMap[dim], eid, index);CHKERRQ(ierr); 135 entities->entity[dim][index].dim = dim; 136 entities->entity[dim][index].id = eid; 137 if (entity) *entity = &entities->entity[dim][index]; 138 PetscFunctionReturn(0); 139 } 140 141 static PetscErrorCode GmshEntitiesGet(GmshEntities *entities, PetscInt dim, PetscInt eid, GmshEntity** entity) 142 { 143 PetscInt index; 144 PetscErrorCode ierr; 145 146 PetscFunctionBegin; 147 ierr = PetscHMapIGet(entities->entityMap[dim], eid, &index);CHKERRQ(ierr); 148 *entity = &entities->entity[dim][index]; 149 PetscFunctionReturn(0); 150 } 151 152 static PetscErrorCode GmshEntitiesDestroy(GmshEntities **entities) 153 { 154 PetscInt dim; 155 PetscErrorCode ierr; 156 157 PetscFunctionBegin; 158 if (!*entities) PetscFunctionReturn(0); 159 for (dim = 0; dim < 4; ++dim) { 160 ierr = PetscFree((*entities)->entity[dim]);CHKERRQ(ierr); 161 ierr = PetscHMapIDestroy(&(*entities)->entityMap[dim]);CHKERRQ(ierr); 162 } 163 ierr = PetscFree((*entities));CHKERRQ(ierr); 164 PetscFunctionReturn(0); 165 } 166 167 typedef struct { 168 PetscInt dim; /* Entity dimension */ 169 PetscInt id; /* Element number */ 170 PetscInt cellType; /* Cell type */ 171 PetscInt numNodes; /* Size of node array */ 172 int nodes[8]; /* Node array */ 173 PetscInt numTags; /* Size of tag array */ 174 int tags[4]; /* Tag array */ 175 } GmshElement; 176 177 static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **coordinates) 178 { 179 char line[PETSC_MAX_PATH_LEN]; 180 int v, nid, snum; 181 PetscErrorCode ierr; 182 183 PetscFunctionBegin; 184 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 185 snum = sscanf(line, "%d", numVertices); 186 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 187 ierr = PetscMalloc1(*numVertices*3, coordinates);CHKERRQ(ierr); 188 for (v = 0; v < *numVertices; ++v) { 189 double *xyz = *coordinates + v*3; 190 ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 191 if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 192 if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 193 ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 194 if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 195 } 196 PetscFunctionReturn(0); 197 } 198 199 /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 200 file contents multiple times to figure out the true number of cells and facets 201 in the given mesh. To make this more efficient we read the file contents only 202 once and store them in memory, while determining the true number of cells. */ 203 static PetscErrorCode DMPlexCreateGmsh_ReadElements_v2(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int *numCells, GmshElement **gmsh_elems) 204 { 205 char line[PETSC_MAX_PATH_LEN]; 206 GmshElement *elements; 207 int i, c, p, ibuf[1+4+512], snum; 208 int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 213 snum = sscanf(line, "%d", numCells); 214 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 215 ierr = PetscMalloc1(*numCells, &elements);CHKERRQ(ierr); 216 for (c = 0; c < *numCells;) { 217 ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 218 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 219 if (binary) { 220 cellType = ibuf[0]; 221 numElem = ibuf[1]; 222 numTags = ibuf[2]; 223 } else { 224 elements[c].id = ibuf[0]; 225 cellType = ibuf[1]; 226 numTags = ibuf[2]; 227 numElem = 1; 228 } 229 /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 230 numNodesIgnore = 0; 231 switch (cellType) { 232 case 1: /* 2-node line */ 233 dim = 1; 234 numNodes = 2; 235 break; 236 case 2: /* 3-node triangle */ 237 dim = 2; 238 numNodes = 3; 239 break; 240 case 3: /* 4-node quadrangle */ 241 dim = 2; 242 numNodes = 4; 243 break; 244 case 4: /* 4-node tetrahedron */ 245 dim = 3; 246 numNodes = 4; 247 break; 248 case 5: /* 8-node hexahedron */ 249 dim = 3; 250 numNodes = 8; 251 break; 252 case 6: /* 6-node wedge */ 253 dim = 3; 254 numNodes = 6; 255 break; 256 case 8: /* 3-node 2nd order line */ 257 dim = 1; 258 numNodes = 2; 259 numNodesIgnore = 1; 260 break; 261 case 9: /* 6-node 2nd order triangle */ 262 dim = 2; 263 numNodes = 3; 264 numNodesIgnore = 3; 265 break; 266 case 13: /* 18-node 2nd wedge */ 267 dim = 3; 268 numNodes = 6; 269 numNodesIgnore = 12; 270 break; 271 case 15: /* 1-node vertex */ 272 dim = 0; 273 numNodes = 1; 274 break; 275 case 7: /* 5-node pyramid */ 276 case 10: /* 9-node 2nd order quadrangle */ 277 case 11: /* 10-node 2nd order tetrahedron */ 278 case 12: /* 27-node 2nd order hexhedron */ 279 case 14: /* 14-node 2nd order pyramid */ 280 default: 281 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 282 } 283 if (binary) { 284 const int nint = 1 + numTags + numNodes + numNodesIgnore; 285 /* Loop over element blocks */ 286 for (i = 0; i < numElem; ++i, ++c) { 287 ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 288 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 289 elements[c].dim = dim; 290 elements[c].numNodes = numNodes; 291 elements[c].numTags = numTags; 292 elements[c].id = ibuf[0]; 293 elements[c].cellType = cellType; 294 for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 295 for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 296 } 297 } else { 298 elements[c].dim = dim; 299 elements[c].numNodes = numNodes; 300 elements[c].numTags = numTags; 301 elements[c].cellType = cellType; 302 ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 303 ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 304 ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 305 c++; 306 } 307 } 308 *gmsh_elems = elements; 309 PetscFunctionReturn(0); 310 } 311 312 /* 313 $Entities 314 numPoints(unsigned long) numCurves(unsigned long) 315 numSurfaces(unsigned long) numVolumes(unsigned long) 316 // points 317 tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 318 boxMaxX(double) boxMaxY(double) boxMaxZ(double) 319 numPhysicals(unsigned long) phyisicalTag[...](int) 320 ... 321 // curves 322 tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 323 boxMaxX(double) boxMaxY(double) boxMaxZ(double) 324 numPhysicals(unsigned long) physicalTag[...](int) 325 numBREPVert(unsigned long) tagBREPVert[...](int) 326 ... 327 // surfaces 328 tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 329 boxMaxX(double) boxMaxY(double) boxMaxZ(double) 330 numPhysicals(unsigned long) physicalTag[...](int) 331 numBREPCurve(unsigned long) tagBREPCurve[...](int) 332 ... 333 // volumes 334 tag(int) boxMinX(double) boxMinY(double) boxMinZ(double) 335 boxMaxX(double) boxMaxY(double) boxMaxZ(double) 336 numPhysicals(unsigned long) physicalTag[...](int) 337 numBREPSurfaces(unsigned long) tagBREPSurfaces[...](int) 338 ... 339 $EndEntities 340 */ 341 static PetscErrorCode DMPlexCreateGmsh_ReadEntities_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshEntities **entities) 342 { 343 long index, num, count[4]; 344 int dim, eid, numTags, *ibuf, t; 345 GmshEntity *entity = NULL; 346 GmshWorkBuffer work; 347 PetscErrorCode ierr; 348 349 PetscFunctionBegin; 350 ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 351 ierr = PetscViewerRead(viewer, count, 4, NULL, PETSC_LONG);CHKERRQ(ierr); 352 if (byteSwap) {ierr = PetscByteSwap(count, PETSC_LONG, 4);CHKERRQ(ierr);} 353 ierr = GmshEntitiesCreate(count, entities);CHKERRQ(ierr); 354 for (dim = 0; dim < 4; ++dim) { 355 for (index = 0; index < count[dim]; ++index) { 356 ierr = PetscViewerRead(viewer, &eid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 357 if (byteSwap) {ierr = PetscByteSwap(&eid, PETSC_ENUM, 1);CHKERRQ(ierr);} 358 ierr = GmshEntitiesAdd(*entities, (PetscInt)index, dim, eid, &entity);CHKERRQ(ierr); 359 ierr = PetscViewerRead(viewer, entity->bbox, 6, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 360 if (byteSwap) {ierr = PetscByteSwap(entity->bbox, PETSC_DOUBLE, 6);CHKERRQ(ierr);} 361 ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 362 if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 363 ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr); 364 ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 365 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 366 entity->numTags = numTags = (int) PetscMin(num, 4); 367 for (t = 0; t < numTags; ++t) entity->tags[t] = ibuf[t]; 368 if (dim == 0) continue; 369 ierr = PetscViewerRead(viewer, &num, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 370 if (byteSwap) {ierr = PetscByteSwap(&num, PETSC_LONG, 1);CHKERRQ(ierr);} 371 ierr = GmshWorkBufferGet(&work, num, sizeof(int), &ibuf);CHKERRQ(ierr); 372 ierr = PetscViewerRead(viewer, ibuf, num, NULL, PETSC_ENUM);CHKERRQ(ierr); 373 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, num);CHKERRQ(ierr);} 374 } 375 } 376 ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 377 PetscFunctionReturn(0); 378 } 379 380 /* 381 $Nodes 382 numEntityBlocks(unsigned long) numNodes(unsigned long) 383 tagEntity(int) dimEntity(int) typeNode(int) numNodes(unsigned long) 384 tag(int) x(double) y(double) z(double) 385 ... 386 ... 387 $EndNodes 388 */ 389 static PetscErrorCode DMPlexCreateGmsh_ReadNodes_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int *numVertices, double **gmsh_nodes) 390 { 391 long block, node, v, numEntityBlocks, numTotalNodes, numNodes; 392 int info[3], nid; 393 double *coordinates; 394 char *cbuf; 395 GmshWorkBuffer work; 396 PetscErrorCode ierr; 397 398 PetscFunctionBegin; 399 ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 400 ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 401 if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 402 ierr = PetscViewerRead(viewer, &numTotalNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 403 if (byteSwap) {ierr = PetscByteSwap(&numTotalNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 404 ierr = PetscMalloc1(numTotalNodes*3, &coordinates);CHKERRQ(ierr); 405 *numVertices = (int)numTotalNodes; 406 *gmsh_nodes = coordinates; 407 for (v = 0, block = 0; block < numEntityBlocks; ++block) { 408 ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 409 ierr = PetscViewerRead(viewer, &numNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 410 if (byteSwap) {ierr = PetscByteSwap(&numNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 411 if (binary) { 412 int nbytes = sizeof(int) + 3*sizeof(double); 413 ierr = GmshWorkBufferGet(&work, numNodes, nbytes, &cbuf);CHKERRQ(ierr); 414 ierr = PetscViewerRead(viewer, cbuf, numNodes*nbytes, NULL, PETSC_CHAR);CHKERRQ(ierr); 415 for (node = 0; node < numNodes; ++node, ++v) { 416 char *cnid = cbuf + node*nbytes, *cxyz = cnid + sizeof(int); 417 double *xyz = coordinates + v*3; 418 #if !defined(PETSC_WORDS_BIGENDIAN) 419 ierr = PetscByteSwap(cnid, PETSC_ENUM, 1);CHKERRQ(ierr); 420 ierr = PetscByteSwap(cxyz, PETSC_DOUBLE, 3);CHKERRQ(ierr); 421 #endif 422 ierr = PetscMemcpy(&nid, cnid, sizeof(int));CHKERRQ(ierr); 423 ierr = PetscMemcpy(xyz, cxyz, 3*sizeof(double));CHKERRQ(ierr); 424 if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 425 if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 426 if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 427 } 428 } else { 429 for (node = 0; node < numNodes; ++node, ++v) { 430 double *xyz = coordinates + v*3; 431 ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 432 if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 433 if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 434 ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 435 if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 436 } 437 } 438 } 439 ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 440 PetscFunctionReturn(0); 441 } 442 443 /* 444 $Elements 445 numEntityBlocks(unsigned long) numElements(unsigned long) 446 tagEntity(int) dimEntity(int) typeEle(int) numElements(unsigned long) 447 tag(int) numVert[...](int) 448 ... 449 ... 450 $EndElements 451 */ 452 static PetscErrorCode DMPlexCreateGmsh_ReadElements_v4(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, GmshEntities *entities, int *numCells, GmshElement **gmsh_elems) 453 { 454 long c, block, numEntityBlocks, numTotalElements, elem, numElements; 455 int p, info[3], *ibuf = NULL; 456 int eid, dim, numTags, *tags, cellType, numNodes; 457 GmshEntity *entity = NULL; 458 GmshElement *elements; 459 GmshWorkBuffer work; 460 PetscErrorCode ierr; 461 462 PetscFunctionBegin; 463 ierr = GmshWorkBufferInit(&work);CHKERRQ(ierr); 464 ierr = PetscViewerRead(viewer, &numEntityBlocks, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 465 if (byteSwap) {ierr = PetscByteSwap(&numEntityBlocks, PETSC_LONG, 1);CHKERRQ(ierr);} 466 ierr = PetscViewerRead(viewer, &numTotalElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 467 if (byteSwap) {ierr = PetscByteSwap(&numTotalElements, PETSC_LONG, 1);CHKERRQ(ierr);} 468 ierr = PetscCalloc1(numTotalElements, &elements);CHKERRQ(ierr); 469 *numCells = (int)numTotalElements; 470 *gmsh_elems = elements; 471 for (c = 0, block = 0; block < numEntityBlocks; ++block) { 472 ierr = PetscViewerRead(viewer, info, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 473 if (byteSwap) {ierr = PetscByteSwap(info, PETSC_ENUM, 3);CHKERRQ(ierr);} 474 eid = info[0]; dim = info[1]; cellType = info[2]; 475 ierr = GmshEntitiesGet(entities, dim, eid, &entity);CHKERRQ(ierr); 476 numTags = entity->numTags; 477 tags = entity->tags; 478 switch (cellType) { 479 case 1: /* 2-node line */ 480 numNodes = 2; 481 break; 482 case 2: /* 3-node triangle */ 483 numNodes = 3; break; 484 case 3: /* 4-node quadrangle */ 485 numNodes = 4; 486 break; 487 case 4: /* 4-node tetrahedron */ 488 numNodes = 4; 489 break; 490 case 5: /* 8-node hexahedron */ 491 numNodes = 8; 492 break; 493 case 6: /* 6-node wedge */ 494 numNodes = 6; 495 break; 496 case 15: /* 1-node vertex */ 497 numNodes = 1; 498 break; 499 default: 500 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 501 } 502 ierr = PetscViewerRead(viewer, &numElements, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 503 if (byteSwap) {ierr = PetscByteSwap(&numElements, PETSC_LONG, 1);CHKERRQ(ierr);} 504 ierr = GmshWorkBufferGet(&work, (1+numNodes)*numElements, sizeof(int), &ibuf);CHKERRQ(ierr); 505 ierr = PetscViewerRead(viewer, ibuf, (1+numNodes)*numElements, NULL, PETSC_ENUM);CHKERRQ(ierr); 506 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, (1+numNodes)*numElements);CHKERRQ(ierr);} 507 for (elem = 0; elem < numElements; ++elem, ++c) { 508 int *id = ibuf + elem*(1+numNodes), *nodes = id + 1; 509 GmshElement *element = elements + c; 510 element->dim = dim; 511 element->cellType = cellType; 512 element->numNodes = numNodes; 513 element->numTags = numTags; 514 element->id = id[0]; 515 for (p = 0; p < numNodes; p++) element->nodes[p] = nodes[p]; 516 for (p = 0; p < numTags; p++) element->tags[p] = tags[p]; 517 } 518 } 519 ierr = GmshWorkBufferFree(&work);CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 523 /*@ 524 DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 525 526 Collective on comm 527 528 Input Parameters: 529 + comm - The MPI communicator 530 . viewer - The Viewer associated with a Gmsh file 531 - interpolate - Create faces and edges in the mesh 532 533 Output Parameter: 534 . dm - The DM object representing the mesh 535 536 Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 537 and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 538 539 Level: beginner 540 541 .keywords: mesh,Gmsh 542 .seealso: DMPLEX, DMCreate() 543 @*/ 544 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 545 { 546 PetscViewer parentviewer = NULL; 547 double *coordsIn = NULL; 548 GmshEntities *entities = NULL; 549 GmshElement *gmsh_elem = NULL; 550 PetscSection coordSection; 551 Vec coordinates; 552 PetscBT periodicV = NULL, periodicC = NULL; 553 PetscScalar *coords; 554 PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 555 int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 556 PetscMPIInt rank; 557 char line[PETSC_MAX_PATH_LEN]; 558 PetscBool binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 559 PetscBool enable_hybrid = PETSC_FALSE; 560 PetscErrorCode ierr; 561 562 PetscFunctionBegin; 563 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 564 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 565 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 566 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 567 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 568 ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 569 if (zerobase) shift = 0; 570 571 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 572 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 573 ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 574 575 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 576 577 /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 578 if (binary) { 579 parentviewer = viewer; 580 ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 581 } 582 583 if (!rank) { 584 PetscBool match, hybrid; 585 int fileType, fileFormat, dataSize; 586 float version; 587 588 /* Read header */ 589 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 590 ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 591 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 592 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 593 snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 594 if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 595 if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 596 if ((int)version == 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file version 3.0 not supported"); 597 if (version > 4.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at most version 4.0"); 598 if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 599 fileFormat = (int)version; 600 if (binary) { 601 int checkInt; 602 ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 603 if (checkInt != 1) { 604 ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 605 if (checkInt == 1) byteSwap = PETSC_TRUE; 606 else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 607 } 608 } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 609 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 610 ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 611 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 612 613 /* OPTIONAL Read physical names */ 614 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 615 ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 616 if (match) { 617 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 618 snum = sscanf(line, "%d", &numRegions); 619 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 620 for (r = 0; r < numRegions; ++r) { 621 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 622 } 623 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 624 ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 625 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 626 /* Initial read for vertex section */ 627 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 628 } 629 630 /* Read entities */ 631 if (fileFormat == 4) { 632 ierr = PetscStrncmp(line, "$Entities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 633 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 634 ierr = DMPlexCreateGmsh_ReadEntities_v4(viewer, binary, byteSwap, &entities);CHKERRQ(ierr); 635 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 636 ierr = PetscStrncmp(line, "$EndEntities", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 637 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 638 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 639 } 640 641 /* Read vertices */ 642 ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 643 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 644 if (fileFormat == 2) { 645 ierr = DMPlexCreateGmsh_ReadNodes_v2(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr); 646 } else { 647 ierr = DMPlexCreateGmsh_ReadNodes_v4(viewer, binary, byteSwap, shift, &numVertices, &coordsIn);CHKERRQ(ierr); 648 } 649 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 650 ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 651 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 652 653 /* Read cells */ 654 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 655 ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 656 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 657 if (fileFormat == 2) { 658 ierr = DMPlexCreateGmsh_ReadElements_v2(viewer, binary, byteSwap, shift, &numCells, &gmsh_elem);CHKERRQ(ierr); 659 } else { 660 ierr = DMPlexCreateGmsh_ReadElements_v4(viewer, binary, byteSwap, shift, entities, &numCells, &gmsh_elem);CHKERRQ(ierr); 661 } 662 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 663 ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 664 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 665 666 ierr = GmshEntitiesDestroy(&entities);CHKERRQ(ierr); 667 668 hybrid = PETSC_FALSE; 669 for (trueNumCells = 0, c = 0; c < numCells; ++c) { 670 int on = -1; 671 if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 672 if (gmsh_elem[c].dim == dim) {hybrid = (trueNumCells ? (on != gmsh_elem[c].cellType ? on = gmsh_elem[c].cellType,PETSC_TRUE : hybrid) : (on = gmsh_elem[c].cellType, PETSC_FALSE) ); trueNumCells++;} 673 /* wedges always indicate an hybrid mesh in PLEX */ 674 if (on == 6 || on == 13) hybrid = PETSC_TRUE; 675 } 676 677 /* Renumber cells for hybrid grids */ 678 if (hybrid && enable_hybrid) { 679 PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 680 PetscInt cell, tn, *tp; 681 int n1 = 0,n2 = 0; 682 683 ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 684 ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 685 for (cell = 0, c = 0; c < numCells; ++c) { 686 if (gmsh_elem[c].dim == dim) { 687 if (!n1) n1 = gmsh_elem[c].cellType; 688 else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 689 690 if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 691 else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 692 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 693 cell++; 694 } 695 } 696 697 switch (n1) { 698 case 2: /* triangles */ 699 case 9: 700 switch (n2) { 701 case 0: /* single type mesh */ 702 case 3: /* quads */ 703 case 10: 704 break; 705 default: 706 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 707 } 708 break; 709 case 3: 710 case 10: 711 switch (n2) { 712 case 0: /* single type mesh */ 713 case 2: /* swap since we list simplices first */ 714 case 9: 715 tn = hc1; 716 hc1 = hc2; 717 hc2 = tn; 718 719 tp = hybridCells1; 720 hybridCells1 = hybridCells2; 721 hybridCells2 = tp; 722 break; 723 default: 724 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 725 } 726 break; 727 case 4: /* tetrahedra */ 728 case 11: 729 switch (n2) { 730 case 0: /* single type mesh */ 731 case 6: /* wedges */ 732 case 13: 733 break; 734 default: 735 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 736 } 737 break; 738 case 6: 739 case 13: 740 switch (n2) { 741 case 0: /* single type mesh */ 742 case 4: /* swap since we list simplices first */ 743 case 11: 744 tn = hc1; 745 hc1 = hc2; 746 hc2 = tn; 747 748 tp = hybridCells1; 749 hybridCells1 = hybridCells2; 750 hybridCells2 = tp; 751 break; 752 default: 753 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 754 } 755 break; 756 default: 757 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 758 } 759 cMax = hc1; 760 ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 761 for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 762 for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 763 ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 764 ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 765 } 766 767 /* OPTIONAL Read periodic section */ 768 if (periodic) { 769 PetscInt pVert, *periodicMapT, *aux; 770 int numPeriodic; 771 772 ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 773 ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 774 for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 775 776 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 777 ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 778 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 779 if (fileFormat == 2 || !binary) { 780 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 781 snum = sscanf(line, "%d", &numPeriodic); 782 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 783 } else { 784 ierr = PetscViewerRead(viewer, &numPeriodic, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 785 if (byteSwap) {ierr = PetscByteSwap(&numPeriodic, PETSC_ENUM, 1);CHKERRQ(ierr);} 786 } 787 for (i = 0; i < numPeriodic; i++) { 788 int ibuf[3], slaveDim = -1, slaveTag = -1, masterTag = -1, slaveNode, masterNode; 789 long j, nNodes; 790 double affine[16]; 791 792 if (fileFormat == 2 || !binary) { 793 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 794 snum = sscanf(line, "%d %d %d", &slaveDim, &slaveTag, &masterTag); 795 if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 796 } else { 797 ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 798 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 799 slaveDim = ibuf[0]; slaveTag = ibuf[1]; masterTag = ibuf[2]; 800 } 801 (void)slaveDim; (void)slaveTag; (void)masterTag; /* unused */ 802 803 if (fileFormat == 2 || !binary) { 804 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 805 snum = sscanf(line, "%ld", &nNodes); 806 if (snum != 1) { /* discard tranformation and try again */ 807 ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 808 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 809 snum = sscanf(line, "%ld", &nNodes); 810 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 811 } 812 } else { 813 ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 814 if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 815 if (nNodes == -1) { /* discard tranformation and try again */ 816 ierr = PetscViewerRead(viewer, affine, 16, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 817 ierr = PetscViewerRead(viewer, &nNodes, 1, NULL, PETSC_LONG);CHKERRQ(ierr); 818 if (byteSwap) {ierr = PetscByteSwap(&nNodes, PETSC_LONG, 1);CHKERRQ(ierr);} 819 } 820 } 821 822 for (j = 0; j < nNodes; j++) { 823 if (fileFormat == 2 || !binary) { 824 ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 825 snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 826 if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 827 } else { 828 ierr = PetscViewerRead(viewer, ibuf, 2, NULL, PETSC_ENUM);CHKERRQ(ierr); 829 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 2);CHKERRQ(ierr);} 830 slaveNode = ibuf[0]; masterNode = ibuf[1]; 831 } 832 periodicMapT[slaveNode - shift] = masterNode - shift; 833 ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr); 834 ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr); 835 } 836 } 837 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 838 ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 839 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 840 841 /* we may have slaves of slaves */ 842 for (i = 0; i < numVertices; i++) { 843 while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 844 periodicMapT[i] = periodicMapT[periodicMapT[i]]; 845 } 846 } 847 /* periodicMap : from old to new numbering (periodic vertices excluded) 848 periodicMapI: from new to old numbering */ 849 ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 850 ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 851 ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 852 for (i = 0, pVert = 0; i < numVertices; i++) { 853 if (periodicMapT[i] != i) { 854 pVert++; 855 } else { 856 aux[i] = i - pVert; 857 periodicMapI[i - pVert] = i; 858 } 859 } 860 for (i = 0 ; i < numVertices; i++) { 861 periodicMap[i] = aux[periodicMapT[i]]; 862 } 863 ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 864 ierr = PetscFree(aux);CHKERRQ(ierr); 865 /* remove periodic vertices */ 866 numVertices = numVertices - pVert; 867 } 868 } 869 870 if (parentviewer) { 871 ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 872 } 873 874 /* Allocate the cell-vertex mesh */ 875 ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 876 for (cell = 0, c = 0; c < numCells; ++c) { 877 if (gmsh_elem[c].dim == dim) { 878 ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 879 cell++; 880 } 881 } 882 ierr = DMSetUp(*dm);CHKERRQ(ierr); 883 /* Add cell-vertex connections */ 884 for (cell = 0, c = 0; c < numCells; ++c) { 885 if (gmsh_elem[c].dim == dim) { 886 PetscInt pcone[8], corner; 887 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 888 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 889 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 890 } 891 if (dim == 3) { 892 /* Tetrahedra are inverted */ 893 if (gmsh_elem[c].cellType == 4) { 894 PetscInt tmp = pcone[0]; 895 pcone[0] = pcone[1]; 896 pcone[1] = tmp; 897 } 898 /* Hexahedra are inverted */ 899 if (gmsh_elem[c].cellType == 5) { 900 PetscInt tmp = pcone[1]; 901 pcone[1] = pcone[3]; 902 pcone[3] = tmp; 903 } 904 /* Prisms are inverted */ 905 if (gmsh_elem[c].cellType == 6) { 906 PetscInt tmp; 907 908 tmp = pcone[1]; 909 pcone[1] = pcone[2]; 910 pcone[2] = tmp; 911 tmp = pcone[4]; 912 pcone[4] = pcone[5]; 913 pcone[5] = tmp; 914 } 915 } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 916 /* quads are input to PLEX as prisms */ 917 if (gmsh_elem[c].cellType == 3) { 918 PetscInt tmp = pcone[2]; 919 pcone[2] = pcone[3]; 920 pcone[3] = tmp; 921 } 922 } 923 ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 924 cell++; 925 } 926 } 927 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 928 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 929 ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 930 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 931 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 932 if (interpolate) { 933 DM idm; 934 935 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 936 ierr = DMDestroy(dm);CHKERRQ(ierr); 937 *dm = idm; 938 } 939 940 if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 941 if (!rank && usemarker) { 942 PetscInt f, fStart, fEnd; 943 944 ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 945 ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 946 for (f = fStart; f < fEnd; ++f) { 947 PetscInt suppSize; 948 949 ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 950 if (suppSize == 1) { 951 PetscInt *cone = NULL, coneSize, p; 952 953 ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 954 for (p = 0; p < coneSize; p += 2) { 955 ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 956 } 957 ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 958 } 959 } 960 } 961 962 if (!rank) { 963 PetscInt vStart, vEnd; 964 965 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 966 for (cell = 0, c = 0; c < numCells; ++c) { 967 968 /* Create face sets */ 969 if (interpolate && gmsh_elem[c].dim == dim-1) { 970 const PetscInt *join; 971 PetscInt joinSize, pcone[8], corner; 972 /* Find the relevant facet with vertex joins */ 973 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 974 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 975 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 976 } 977 ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 978 if (joinSize != 1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for GMsh element %d (Plex cell %D)", gmsh_elem[c].id, c); 979 ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 980 ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 981 } 982 983 /* Create cell sets */ 984 if (gmsh_elem[c].dim == dim) { 985 if (gmsh_elem[c].numTags > 0) { 986 ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 987 } 988 cell++; 989 } 990 991 /* Create vertex sets */ 992 if (gmsh_elem[c].dim == 0) { 993 if (gmsh_elem[c].numTags > 0) { 994 const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 995 const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 996 ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 997 } 998 } 999 } 1000 } 1001 1002 /* Create coordinates */ 1003 if (embedDim < 0) embedDim = dim; 1004 ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 1005 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1006 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1007 ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 1008 if (periodicMap) { /* we need to localize coordinates on cells */ 1009 ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 1010 } else { 1011 ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 1012 } 1013 for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 1014 ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 1015 ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 1016 } 1017 if (periodicMap) { 1018 ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 1019 for (cell = 0, c = 0; c < numCells; ++c) { 1020 if (gmsh_elem[c].dim == dim) { 1021 PetscInt corner; 1022 PetscBool pc = PETSC_FALSE; 1023 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1024 pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 1025 } 1026 if (pc) { 1027 PetscInt dof = gmsh_elem[c].numNodes*embedDim; 1028 PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1029 ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 1030 ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 1031 ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 1032 } 1033 cell++; 1034 } 1035 } 1036 } 1037 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1038 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 1039 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1040 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1041 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 1042 ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 1043 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 1044 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 1045 if (periodicMap) { 1046 PetscInt off; 1047 1048 for (cell = 0, c = 0; c < numCells; ++c) { 1049 PetscInt pcone[8], corner; 1050 if (gmsh_elem[c].dim == dim) { 1051 PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 1052 if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 1053 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1054 pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 1055 } 1056 if (dim == 3) { 1057 /* Tetrahedra are inverted */ 1058 if (gmsh_elem[c].cellType == 4) { 1059 PetscInt tmp = pcone[0]; 1060 pcone[0] = pcone[1]; 1061 pcone[1] = tmp; 1062 } 1063 /* Hexahedra are inverted */ 1064 if (gmsh_elem[c].cellType == 5) { 1065 PetscInt tmp = pcone[1]; 1066 pcone[1] = pcone[3]; 1067 pcone[3] = tmp; 1068 } 1069 /* Prisms are inverted */ 1070 if (gmsh_elem[c].cellType == 6) { 1071 PetscInt tmp; 1072 1073 tmp = pcone[1]; 1074 pcone[1] = pcone[2]; 1075 pcone[2] = tmp; 1076 tmp = pcone[4]; 1077 pcone[4] = pcone[5]; 1078 pcone[5] = tmp; 1079 } 1080 } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 1081 /* quads are input to PLEX as prisms */ 1082 if (gmsh_elem[c].cellType == 3) { 1083 PetscInt tmp = pcone[2]; 1084 pcone[2] = pcone[3]; 1085 pcone[3] = tmp; 1086 } 1087 } 1088 ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 1089 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 1090 v = pcone[corner]; 1091 for (d = 0; d < embedDim; ++d) { 1092 coords[off++] = (PetscReal) coordsIn[v*3+d]; 1093 } 1094 } 1095 } 1096 cell++; 1097 } 1098 } 1099 for (v = 0; v < numVertices; ++v) { 1100 ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 1101 for (d = 0; d < embedDim; ++d) { 1102 coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 1103 } 1104 } 1105 } else { 1106 for (v = 0; v < numVertices; ++v) { 1107 for (d = 0; d < embedDim; ++d) { 1108 coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 1109 } 1110 } 1111 } 1112 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1113 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1114 ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 1115 1116 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 1117 ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 1118 ierr = PetscFree(hybridMap);CHKERRQ(ierr); 1119 ierr = PetscFree(periodicMap);CHKERRQ(ierr); 1120 ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 1121 ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 1122 ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 1123 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1124 1125 ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 1126 PetscFunctionReturn(0); 1127 } 1128