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