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 ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 50 } 51 ierr = MPI_Bcast(&fileType, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 52 vtype = (fileType == 0) ? PETSCVIEWERASCII : PETSCVIEWERBINARY; 53 54 /* Create appropriate viewer and build plex */ 55 ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 56 ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 57 ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 58 ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 59 ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 60 ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 61 PetscFunctionReturn(0); 62 } 63 64 static PetscErrorCode DMPlexCreateGmsh_ReadNodes(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, int shift, int numVertices, double **coordinates) 65 { 66 int v,nid; 67 PetscErrorCode ierr; 68 69 PetscFunctionBegin; 70 ierr = PetscMalloc1(numVertices*3, coordinates);CHKERRQ(ierr); 71 for (v = 0; v < numVertices; ++v) { 72 double *xyz = *coordinates + v*3; 73 ierr = PetscViewerRead(viewer, &nid, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 74 if (byteSwap) {ierr = PetscByteSwap(&nid, PETSC_ENUM, 1);CHKERRQ(ierr);} 75 if (nid != v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unexpected node number %d should be %d", nid, v+shift); 76 ierr = PetscViewerRead(viewer, xyz, 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 77 if (byteSwap) {ierr = PetscByteSwap(xyz, PETSC_DOUBLE, 3);CHKERRQ(ierr);} 78 } 79 PetscFunctionReturn(0); 80 } 81 82 static PetscErrorCode DMPlexCreateGmsh_ReadElements(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, PETSC_UNUSED int shift, int numCells, GmshElement **gmsh_elems) 83 { 84 GmshElement *elements; 85 int i, c, p, ibuf[1+4+512]; 86 int cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 87 PetscErrorCode ierr; 88 89 PetscFunctionBegin; 90 ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 91 for (c = 0; c < numCells;) { 92 ierr = PetscViewerRead(viewer, ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 93 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, 3);CHKERRQ(ierr);} 94 if (binary) { 95 cellType = ibuf[0]; 96 numElem = ibuf[1]; 97 numTags = ibuf[2]; 98 } else { 99 elements[c].id = ibuf[0]; 100 cellType = ibuf[1]; 101 numTags = ibuf[2]; 102 numElem = 1; 103 } 104 /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 105 numNodesIgnore = 0; 106 switch (cellType) { 107 case 1: /* 2-node line */ 108 dim = 1; 109 numNodes = 2; 110 break; 111 case 2: /* 3-node triangle */ 112 dim = 2; 113 numNodes = 3; 114 break; 115 case 3: /* 4-node quadrangle */ 116 dim = 2; 117 numNodes = 4; 118 break; 119 case 4: /* 4-node tetrahedron */ 120 dim = 3; 121 numNodes = 4; 122 break; 123 case 5: /* 8-node hexahedron */ 124 dim = 3; 125 numNodes = 8; 126 break; 127 case 6: /* 6-node wedge */ 128 dim = 3; 129 numNodes = 6; 130 break; 131 case 8: /* 3-node 2nd order line */ 132 dim = 1; 133 numNodes = 2; 134 numNodesIgnore = 1; 135 break; 136 case 9: /* 6-node 2nd order triangle */ 137 dim = 2; 138 numNodes = 3; 139 numNodesIgnore = 3; 140 break; 141 case 13: /* 18-node 2nd wedge */ 142 dim = 3; 143 numNodes = 6; 144 numNodesIgnore = 12; 145 break; 146 case 15: /* 1-node vertex */ 147 dim = 0; 148 numNodes = 1; 149 break; 150 case 7: /* 5-node pyramid */ 151 case 10: /* 9-node 2nd order quadrangle */ 152 case 11: /* 10-node 2nd order tetrahedron */ 153 case 12: /* 27-node 2nd order hexhedron */ 154 case 14: /* 14-node 2nd order pyramid */ 155 default: 156 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 157 } 158 if (binary) { 159 const int nint = 1 + numTags + numNodes + numNodesIgnore; 160 /* Loop over element blocks */ 161 for (i = 0; i < numElem; ++i, ++c) { 162 ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 163 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 164 elements[c].dim = dim; 165 elements[c].numNodes = numNodes; 166 elements[c].numTags = numTags; 167 elements[c].id = ibuf[0]; 168 elements[c].cellType = cellType; 169 for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 170 for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 171 } 172 } else { 173 elements[c].dim = dim; 174 elements[c].numNodes = numNodes; 175 elements[c].numTags = numTags; 176 elements[c].cellType = cellType; 177 ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 178 ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 179 ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 180 c++; 181 } 182 } 183 *gmsh_elems = elements; 184 PetscFunctionReturn(0); 185 } 186 187 /*@ 188 DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 189 190 Collective on comm 191 192 Input Parameters: 193 + comm - The MPI communicator 194 . viewer - The Viewer associated with a Gmsh file 195 - interpolate - Create faces and edges in the mesh 196 197 Output Parameter: 198 . dm - The DM object representing the mesh 199 200 Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 201 and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 202 203 Level: beginner 204 205 .keywords: mesh,Gmsh 206 .seealso: DMPLEX, DMCreate() 207 @*/ 208 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 209 { 210 PetscViewer parentviewer = NULL; 211 double *coordsIn = NULL; 212 GmshElement *gmsh_elem = NULL; 213 PetscSection coordSection; 214 Vec coordinates; 215 PetscBT periodicV = NULL, periodicC = NULL; 216 PetscScalar *coords; 217 PetscInt dim = 0, embedDim = -1, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL, *hybridMap = NULL, cMax = PETSC_DETERMINE; 218 int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 219 PetscMPIInt rank; 220 char line[PETSC_MAX_PATH_LEN]; 221 PetscBool binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 222 PetscBool enable_hybrid = PETSC_FALSE; 223 PetscErrorCode ierr; 224 225 PetscFunctionBegin; 226 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 227 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_hybrid", &enable_hybrid, NULL);CHKERRQ(ierr); 228 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 229 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 230 ierr = PetscOptionsGetBool(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 231 ierr = PetscOptionsGetInt(((PetscObject) viewer)->options,((PetscObject) viewer)->prefix, "-dm_plex_gmsh_spacedim", &embedDim, NULL);CHKERRQ(ierr); 232 if (zerobase) shift = 0; 233 234 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 235 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 236 ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 237 238 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 239 240 /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 241 if (binary) { 242 parentviewer = viewer; 243 ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 244 } 245 246 if (!rank) { 247 PetscBool match, hybrid; 248 int fileType, dataSize; 249 float version; 250 251 /* Read header */ 252 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 253 ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 254 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 255 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 256 snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 257 if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 258 if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 259 if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 260 if (binary) { 261 int checkInt; 262 ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 263 if (checkInt != 1) { 264 ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 265 if (checkInt == 1) byteSwap = PETSC_TRUE; 266 else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 267 } 268 } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 269 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 270 ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 271 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 272 273 /* OPTIONAL Read physical names */ 274 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 275 ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 276 if (match) { 277 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 278 snum = sscanf(line, "%d", &numRegions); 279 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 280 for (r = 0; r < numRegions; ++r) { 281 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 282 } 283 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 284 ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 285 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 286 /* Initial read for vertex section */ 287 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 288 } 289 290 /* Read vertices */ 291 ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 292 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 293 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 294 snum = sscanf(line, "%d", &numVertices); 295 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 296 ierr = DMPlexCreateGmsh_ReadNodes(viewer, binary, byteSwap, shift, numVertices, &coordsIn);CHKERRQ(ierr); 297 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 298 ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 299 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 300 301 /* Read cells */ 302 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 303 ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 304 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 305 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 306 snum = sscanf(line, "%d", &numCells); 307 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 308 /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 309 file contents multiple times to figure out the true number of cells and facets 310 in the given mesh. To make this more efficient we read the file contents only 311 once and store them in memory, while determining the true number of cells. */ 312 ierr = DMPlexCreateGmsh_ReadElements(viewer, binary, byteSwap, shift, numCells, &gmsh_elem);CHKERRQ(ierr); 313 hybrid = PETSC_FALSE; 314 for (trueNumCells = 0, c = 0; c < numCells; ++c) { 315 int on = -1; 316 if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 317 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++;} 318 /* wedges always indicate an hybrid mesh in PLEX */ 319 if (on == 6 || on == 13) hybrid = PETSC_TRUE; 320 } 321 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 322 ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 323 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 324 325 /* Renumber cells for hybrid grids */ 326 if (hybrid && enable_hybrid) { 327 PetscInt hc1 = 0, hc2 = 0, *hybridCells1 = NULL, *hybridCells2 = NULL; 328 PetscInt cell, tn, *tp; 329 int n1 = 0,n2 = 0; 330 331 ierr = PetscMalloc1(trueNumCells, &hybridCells1);CHKERRQ(ierr); 332 ierr = PetscMalloc1(trueNumCells, &hybridCells2);CHKERRQ(ierr); 333 for (cell = 0, c = 0; c < numCells; ++c) { 334 if (gmsh_elem[c].dim == dim) { 335 if (!n1) n1 = gmsh_elem[c].cellType; 336 else if (!n2 && n1 != gmsh_elem[c].cellType) n2 = gmsh_elem[c].cellType; 337 338 if (gmsh_elem[c].cellType == n1) hybridCells1[hc1++] = cell; 339 else if (gmsh_elem[c].cellType == n2) hybridCells2[hc2++] = cell; 340 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle more than 2 cell types"); 341 cell++; 342 } 343 } 344 345 switch (n1) { 346 case 2: /* triangles */ 347 case 9: 348 switch (n2) { 349 case 0: /* single type mesh */ 350 case 3: /* quads */ 351 case 10: 352 break; 353 default: 354 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 355 } 356 break; 357 case 3: 358 case 10: 359 switch (n2) { 360 case 0: /* single type mesh */ 361 case 2: /* swap since we list simplices first */ 362 case 9: 363 tn = hc1; 364 hc1 = hc2; 365 hc2 = tn; 366 367 tp = hybridCells1; 368 hybridCells1 = hybridCells2; 369 hybridCells2 = tp; 370 break; 371 default: 372 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 373 } 374 break; 375 case 4: /* tetrahedra */ 376 case 11: 377 switch (n2) { 378 case 0: /* single type mesh */ 379 case 6: /* wedges */ 380 case 13: 381 break; 382 default: 383 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 384 } 385 break; 386 case 6: 387 case 13: 388 switch (n2) { 389 case 0: /* single type mesh */ 390 case 4: /* swap since we list simplices first */ 391 case 11: 392 tn = hc1; 393 hc1 = hc2; 394 hc2 = tn; 395 396 tp = hybridCells1; 397 hybridCells1 = hybridCells2; 398 hybridCells2 = tp; 399 break; 400 default: 401 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 402 } 403 break; 404 default: 405 SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_SUP, "Unsupported cell types %d and %d",n1, n2); 406 } 407 cMax = hc1; 408 ierr = PetscMalloc1(trueNumCells, &hybridMap);CHKERRQ(ierr); 409 for (cell = 0; cell < hc1; cell++) hybridMap[hybridCells1[cell]] = cell; 410 for (cell = 0; cell < hc2; cell++) hybridMap[hybridCells2[cell]] = cell + hc1; 411 ierr = PetscFree(hybridCells1);CHKERRQ(ierr); 412 ierr = PetscFree(hybridCells2);CHKERRQ(ierr); 413 } 414 415 /* OPTIONAL Read periodic section */ 416 if (periodic) { 417 PetscInt pVert, *periodicMapT, *aux; 418 int numPeriodic; 419 420 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 421 ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 422 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 423 ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 424 ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 425 for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 426 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 427 snum = sscanf(line, "%d", &numPeriodic); 428 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 429 for (i = 0; i < numPeriodic; i++) { 430 int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes, slaveNode, masterNode; 431 432 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 433 snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */ 434 if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 435 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 436 snum = sscanf(line, "%d", &nNodes); 437 if (snum != 1) { /* discard tranformation and try again */ 438 ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 439 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 440 snum = sscanf(line, "%d", &nNodes); 441 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 442 } 443 for (j = 0; j < nNodes; j++) { 444 ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 445 snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 446 if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 447 periodicMapT[slaveNode - shift] = masterNode - shift; 448 ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr); 449 ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr); 450 } 451 } 452 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 453 ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 454 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 455 /* we may have slaves of slaves */ 456 for (i = 0; i < numVertices; i++) { 457 while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 458 periodicMapT[i] = periodicMapT[periodicMapT[i]]; 459 } 460 } 461 /* periodicMap : from old to new numbering (periodic vertices excluded) 462 periodicMapI: from new to old numbering */ 463 ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 464 ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 465 ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 466 for (i = 0, pVert = 0; i < numVertices; i++) { 467 if (periodicMapT[i] != i) { 468 pVert++; 469 } else { 470 aux[i] = i - pVert; 471 periodicMapI[i - pVert] = i; 472 } 473 } 474 for (i = 0 ; i < numVertices; i++) { 475 periodicMap[i] = aux[periodicMapT[i]]; 476 } 477 ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 478 ierr = PetscFree(aux);CHKERRQ(ierr); 479 /* remove periodic vertices */ 480 numVertices = numVertices - pVert; 481 } 482 } 483 484 if (parentviewer) { 485 ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 486 } 487 488 /* Allocate the cell-vertex mesh */ 489 ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 490 for (cell = 0, c = 0; c < numCells; ++c) { 491 if (gmsh_elem[c].dim == dim) { 492 ierr = DMPlexSetConeSize(*dm, hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 493 cell++; 494 } 495 } 496 ierr = DMSetUp(*dm);CHKERRQ(ierr); 497 /* Add cell-vertex connections */ 498 for (cell = 0, c = 0; c < numCells; ++c) { 499 if (gmsh_elem[c].dim == dim) { 500 PetscInt pcone[8], corner; 501 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 502 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 503 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 504 } 505 if (dim == 3) { 506 /* Tetrahedra are inverted */ 507 if (gmsh_elem[c].cellType == 4) { 508 PetscInt tmp = pcone[0]; 509 pcone[0] = pcone[1]; 510 pcone[1] = tmp; 511 } 512 /* Hexahedra are inverted */ 513 if (gmsh_elem[c].cellType == 5) { 514 PetscInt tmp = pcone[1]; 515 pcone[1] = pcone[3]; 516 pcone[3] = tmp; 517 } 518 /* Prisms are inverted */ 519 if (gmsh_elem[c].cellType == 6) { 520 PetscInt tmp; 521 522 tmp = pcone[1]; 523 pcone[1] = pcone[2]; 524 pcone[2] = tmp; 525 tmp = pcone[4]; 526 pcone[4] = pcone[5]; 527 pcone[5] = tmp; 528 } 529 } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 530 /* quads are input to PLEX as prisms */ 531 if (gmsh_elem[c].cellType == 3) { 532 PetscInt tmp = pcone[2]; 533 pcone[2] = pcone[3]; 534 pcone[3] = tmp; 535 } 536 } 537 ierr = DMPlexSetCone(*dm, hybridMap ? hybridMap[cell] : cell, pcone);CHKERRQ(ierr); 538 cell++; 539 } 540 } 541 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 542 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 543 ierr = DMPlexSetHybridBounds(*dm, cMax, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 544 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 545 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 546 if (interpolate) { 547 DM idm; 548 549 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 550 ierr = DMDestroy(dm);CHKERRQ(ierr); 551 *dm = idm; 552 } 553 554 if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 555 if (!rank && usemarker) { 556 PetscInt f, fStart, fEnd; 557 558 ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 559 ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 560 for (f = fStart; f < fEnd; ++f) { 561 PetscInt suppSize; 562 563 ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 564 if (suppSize == 1) { 565 PetscInt *cone = NULL, coneSize, p; 566 567 ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 568 for (p = 0; p < coneSize; p += 2) { 569 ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 570 } 571 ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 572 } 573 } 574 } 575 576 if (!rank) { 577 PetscInt vStart, vEnd; 578 579 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 580 for (cell = 0, c = 0; c < numCells; ++c) { 581 582 /* Create face sets */ 583 if (interpolate && gmsh_elem[c].dim == dim-1) { 584 const PetscInt *join; 585 PetscInt joinSize, pcone[8], corner; 586 /* Find the relevant facet with vertex joins */ 587 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 588 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 589 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 590 } 591 ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 592 if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 593 ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 594 ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 595 } 596 597 /* Create cell sets */ 598 if (gmsh_elem[c].dim == dim) { 599 if (gmsh_elem[c].numTags > 0) { 600 ierr = DMSetLabelValue(*dm, "Cell Sets", hybridMap ? hybridMap[cell] : cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 601 } 602 cell++; 603 } 604 605 /* Create vertex sets */ 606 if (gmsh_elem[c].dim == 0) { 607 if (gmsh_elem[c].numTags > 0) { 608 const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 609 const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 610 ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 611 } 612 } 613 } 614 } 615 616 /* Create coordinates */ 617 if (embedDim < 0) embedDim = dim; 618 ierr = DMSetCoordinateDim(*dm, embedDim);CHKERRQ(ierr); 619 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 620 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 621 ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 622 if (periodicMap) { /* we need to localize coordinates on cells */ 623 ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 624 } else { 625 ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 626 } 627 for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 628 ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 629 ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 630 } 631 if (periodicMap) { 632 ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 633 for (cell = 0, c = 0; c < numCells; ++c) { 634 if (gmsh_elem[c].dim == dim) { 635 PetscInt corner; 636 PetscBool pc = PETSC_FALSE; 637 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 638 pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 639 } 640 if (pc) { 641 PetscInt dof = gmsh_elem[c].numNodes*embedDim; 642 PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 643 ierr = PetscBTSet(periodicC, ucell);CHKERRQ(ierr); 644 ierr = PetscSectionSetDof(coordSection, ucell, dof);CHKERRQ(ierr); 645 ierr = PetscSectionSetFieldDof(coordSection, ucell, 0, dof);CHKERRQ(ierr); 646 } 647 cell++; 648 } 649 } 650 } 651 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 652 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 653 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 654 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 655 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 656 ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 657 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 658 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 659 if (periodicMap) { 660 PetscInt off; 661 662 for (cell = 0, c = 0; c < numCells; ++c) { 663 PetscInt pcone[8], corner; 664 if (gmsh_elem[c].dim == dim) { 665 PetscInt ucell = hybridMap ? hybridMap[cell] : cell; 666 if (PetscUnlikely(PetscBTLookup(periodicC, ucell))) { 667 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 668 pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 669 } 670 if (dim == 3) { 671 /* Tetrahedra are inverted */ 672 if (gmsh_elem[c].cellType == 4) { 673 PetscInt tmp = pcone[0]; 674 pcone[0] = pcone[1]; 675 pcone[1] = tmp; 676 } 677 /* Hexahedra are inverted */ 678 if (gmsh_elem[c].cellType == 5) { 679 PetscInt tmp = pcone[1]; 680 pcone[1] = pcone[3]; 681 pcone[3] = tmp; 682 } 683 /* Prisms are inverted */ 684 if (gmsh_elem[c].cellType == 6) { 685 PetscInt tmp; 686 687 tmp = pcone[1]; 688 pcone[1] = pcone[2]; 689 pcone[2] = tmp; 690 tmp = pcone[4]; 691 pcone[4] = pcone[5]; 692 pcone[5] = tmp; 693 } 694 } else if (dim == 2 && hybridMap && hybridMap[cell] >= cMax) { /* hybrid cells */ 695 /* quads are input to PLEX as prisms */ 696 if (gmsh_elem[c].cellType == 3) { 697 PetscInt tmp = pcone[2]; 698 pcone[2] = pcone[3]; 699 pcone[3] = tmp; 700 } 701 } 702 ierr = PetscSectionGetOffset(coordSection, ucell, &off);CHKERRQ(ierr); 703 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 704 v = pcone[corner]; 705 for (d = 0; d < embedDim; ++d) { 706 coords[off++] = (PetscReal) coordsIn[v*3+d]; 707 } 708 } 709 } 710 cell++; 711 } 712 } 713 for (v = 0; v < numVertices; ++v) { 714 ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 715 for (d = 0; d < embedDim; ++d) { 716 coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 717 } 718 } 719 } else { 720 for (v = 0; v < numVertices; ++v) { 721 for (d = 0; d < embedDim; ++d) { 722 coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 723 } 724 } 725 } 726 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 727 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 728 ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 729 730 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 731 ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 732 ierr = PetscFree(hybridMap);CHKERRQ(ierr); 733 ierr = PetscFree(periodicMap);CHKERRQ(ierr); 734 ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 735 ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 736 ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 737 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 738 739 ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 740 PetscFunctionReturn(0); 741 } 742