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 8: /* 3-node 2nd order line */ 128 dim = 1; 129 numNodes = 2; 130 numNodesIgnore = 1; 131 break; 132 case 9: /* 6-node 2nd order triangle */ 133 dim = 2; 134 numNodes = 3; 135 numNodesIgnore = 3; 136 break; 137 case 15: /* 1-node vertex */ 138 dim = 0; 139 numNodes = 1; 140 break; 141 case 6: /* 6-node prism */ 142 case 7: /* 5-node pyramid */ 143 case 10: /* 9-node 2nd order quadrangle */ 144 case 11: /* 10-node 2nd order tetrahedron */ 145 case 12: /* 27-node 2nd order hexhedron */ 146 case 13: /* 19-node 2nd order prism */ 147 case 14: /* 14-node 2nd order pyramid */ 148 default: 149 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 150 } 151 if (binary) { 152 const int nint = 1 + numTags + numNodes + numNodesIgnore; 153 /* Loop over element blocks */ 154 for (i = 0; i < numElem; ++i, ++c) { 155 ierr = PetscViewerRead(viewer, ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 156 if (byteSwap) {ierr = PetscByteSwap(ibuf, PETSC_ENUM, nint);CHKERRQ(ierr);} 157 elements[c].dim = dim; 158 elements[c].numNodes = numNodes; 159 elements[c].numTags = numTags; 160 elements[c].id = ibuf[0]; 161 for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 162 for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 163 } 164 } else { 165 elements[c].dim = dim; 166 elements[c].numNodes = numNodes; 167 elements[c].numTags = numTags; 168 ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 169 ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 170 ierr = PetscViewerRead(viewer, ibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 171 c++; 172 } 173 } 174 *gmsh_elems = elements; 175 PetscFunctionReturn(0); 176 } 177 178 /*@ 179 DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 180 181 Collective on comm 182 183 Input Parameters: 184 + comm - The MPI communicator 185 . viewer - The Viewer associated with a Gmsh file 186 - interpolate - Create faces and edges in the mesh 187 188 Output Parameter: 189 . dm - The DM object representing the mesh 190 191 Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 192 and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 193 194 Level: beginner 195 196 .keywords: mesh,Gmsh 197 .seealso: DMPLEX, DMCreate() 198 @*/ 199 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 200 { 201 PetscViewer parentviewer = NULL; 202 double *coordsIn = NULL; 203 GmshElement *gmsh_elem = NULL; 204 PetscSection coordSection; 205 Vec coordinates; 206 PetscBT periodicV = NULL, periodicC = NULL; 207 PetscScalar *coords; 208 PetscInt dim = 0, embedDim, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL; 209 int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 210 PetscMPIInt rank; 211 char line[PETSC_MAX_PATH_LEN]; 212 PetscBool binary, byteSwap = PETSC_FALSE, zerobase = PETSC_FALSE, isbd = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 213 PetscErrorCode ierr; 214 215 PetscFunctionBegin; 216 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 217 ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_gmsh_periodic", &periodic, NULL);CHKERRQ(ierr); 218 ierr = PetscOptionsGetBool(NULL, NULL, "-dm_plex_gmsh_use_marker", &usemarker, NULL);CHKERRQ(ierr); 219 ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 220 ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr); 221 if (zerobase) shift = 0; 222 223 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 224 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 225 ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 226 227 ierr = PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 228 229 /* Binary viewers read on all ranks, get subviewer to read only in rank 0 */ 230 if (binary) { 231 parentviewer = viewer; 232 ierr = PetscViewerGetSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 233 } 234 235 if (!rank) { 236 PetscBool match; 237 int fileType, dataSize; 238 float version; 239 240 /* Read header */ 241 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 242 ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 243 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 244 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 245 snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 246 if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 247 if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 248 if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 249 if (binary) { 250 int checkInt; 251 ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 252 if (checkInt != 1) { 253 ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 254 if (checkInt == 1) byteSwap = PETSC_TRUE; 255 else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 256 } 257 } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 258 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 259 ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 260 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 261 262 /* OPTIONAL Read physical names */ 263 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 264 ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 265 if (match) { 266 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 267 snum = sscanf(line, "%d", &numRegions); 268 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 269 for (r = 0; r < numRegions; ++r) { 270 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 271 } 272 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 273 ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 274 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 275 /* Initial read for vertex section */ 276 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 277 } 278 279 /* Read vertices */ 280 ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 281 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 282 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 283 snum = sscanf(line, "%d", &numVertices); 284 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 285 ierr = DMPlexCreateGmsh_ReadNodes(viewer, binary, byteSwap, shift, numVertices, &coordsIn);CHKERRQ(ierr); 286 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 287 ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 288 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 289 290 /* Read cells */ 291 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 292 ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 293 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 294 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 295 snum = sscanf(line, "%d", &numCells); 296 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 297 /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 298 file contents multiple times to figure out the true number of cells and facets 299 in the given mesh. To make this more efficient we read the file contents only 300 once and store them in memory, while determining the true number of cells. */ 301 ierr = DMPlexCreateGmsh_ReadElements(viewer, binary, byteSwap, shift, numCells, &gmsh_elem);CHKERRQ(ierr); 302 for (trueNumCells=0, c = 0; c < numCells; ++c) { 303 if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 304 if (gmsh_elem[c].dim == dim) trueNumCells++; 305 } 306 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 307 ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 308 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 309 310 /* OPTIONAL Read periodic section */ 311 if (periodic) { 312 PetscInt pVert, *periodicMapT, *aux; 313 int numPeriodic; 314 315 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 316 ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 317 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 318 ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 319 ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 320 for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 321 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 322 snum = sscanf(line, "%d", &numPeriodic); 323 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 324 for (i = 0; i < numPeriodic; i++) { 325 int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes, slaveNode, masterNode; 326 327 ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 328 snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */ 329 if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 330 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 331 snum = sscanf(line, "%d", &nNodes); 332 if (snum != 1) { /* discard tranformation and try again */ 333 ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 334 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 335 snum = sscanf(line, "%d", &nNodes); 336 if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 337 } 338 for (j = 0; j < nNodes; j++) { 339 ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 340 snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 341 if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 342 periodicMapT[slaveNode - shift] = masterNode - shift; 343 ierr = PetscBTSet(periodicV, slaveNode - shift);CHKERRQ(ierr); 344 ierr = PetscBTSet(periodicV, masterNode - shift);CHKERRQ(ierr); 345 } 346 } 347 ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 348 ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 349 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 350 /* we may have slaves of slaves */ 351 for (i = 0; i < numVertices; i++) { 352 while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 353 periodicMapT[i] = periodicMapT[periodicMapT[i]]; 354 } 355 } 356 /* periodicMap : from old to new numbering (periodic vertices excluded) 357 periodicMapI: from new to old numbering */ 358 ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 359 ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 360 ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 361 for (i = 0, pVert = 0; i < numVertices; i++) { 362 if (periodicMapT[i] != i) { 363 pVert++; 364 } else { 365 aux[i] = i - pVert; 366 periodicMapI[i - pVert] = i; 367 } 368 } 369 for (i = 0 ; i < numVertices; i++) { 370 periodicMap[i] = aux[periodicMapT[i]]; 371 } 372 ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 373 ierr = PetscFree(aux);CHKERRQ(ierr); 374 /* remove periodic vertices */ 375 numVertices = numVertices - pVert; 376 } 377 } 378 379 if (parentviewer) { 380 ierr = PetscViewerRestoreSubViewer(parentviewer, PETSC_COMM_SELF, &viewer);CHKERRQ(ierr); 381 } 382 383 /* Allocate the cell-vertex mesh */ 384 ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 385 for (cell = 0, c = 0; c < numCells; ++c) { 386 if (gmsh_elem[c].dim == dim) { 387 ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 388 cell++; 389 } 390 } 391 ierr = DMSetUp(*dm);CHKERRQ(ierr); 392 /* Add cell-vertex connections */ 393 for (cell = 0, c = 0; c < numCells; ++c) { 394 if (gmsh_elem[c].dim == dim) { 395 PetscInt pcone[8], corner; 396 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 397 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 398 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + trueNumCells; 399 } 400 if (dim == 3) { 401 /* Tetrahedra are inverted */ 402 if (gmsh_elem[c].numNodes == 4) { 403 PetscInt tmp = pcone[0]; 404 pcone[0] = pcone[1]; 405 pcone[1] = tmp; 406 } 407 /* Hexahedra are inverted */ 408 if (gmsh_elem[c].numNodes == 8) { 409 PetscInt tmp = pcone[1]; 410 pcone[1] = pcone[3]; 411 pcone[3] = tmp; 412 } 413 } 414 ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 415 cell++; 416 } 417 } 418 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 419 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 420 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 421 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 422 if (interpolate) { 423 DM idm; 424 425 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 426 ierr = DMDestroy(dm);CHKERRQ(ierr); 427 *dm = idm; 428 } 429 430 if (usemarker && !interpolate && dim > 1) SETERRQ(comm,PETSC_ERR_SUP,"Cannot create marker label without interpolation"); 431 if (!rank && usemarker) { 432 PetscInt f, fStart, fEnd; 433 434 ierr = DMCreateLabel(*dm, "marker");CHKERRQ(ierr); 435 ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 436 for (f = fStart; f < fEnd; ++f) { 437 PetscInt suppSize; 438 439 ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 440 if (suppSize == 1) { 441 PetscInt *cone = NULL, coneSize, p; 442 443 ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 444 for (p = 0; p < coneSize; p += 2) { 445 ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 446 } 447 ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 448 } 449 } 450 } 451 452 if (!rank) { 453 PetscInt vStart, vEnd; 454 455 ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 456 for (cell = 0, c = 0; c < numCells; ++c) { 457 458 /* Create face sets */ 459 if (interpolate && gmsh_elem[c].dim == dim-1) { 460 const PetscInt *join; 461 PetscInt joinSize, pcone[8], corner; 462 /* Find the relevant facet with vertex joins */ 463 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 464 const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 465 pcone[corner] = (periodicMap ? periodicMap[cc] : cc) + vStart; 466 } 467 ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, pcone, &joinSize, &join);CHKERRQ(ierr); 468 if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 469 ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 470 ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 471 } 472 473 /* Create cell sets */ 474 if (gmsh_elem[c].dim == dim) { 475 if (gmsh_elem[c].numTags > 0) { 476 ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 477 } 478 cell++; 479 } 480 481 /* Create vertex sets */ 482 if (gmsh_elem[c].dim == 0) { 483 if (gmsh_elem[c].numTags > 0) { 484 const PetscInt cc = gmsh_elem[c].nodes[0] - shift; 485 const PetscInt vid = (periodicMap ? periodicMap[cc] : cc) + vStart; 486 ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 487 } 488 } 489 } 490 } 491 492 /* Create coordinates */ 493 embedDim = isbd ? dim+1 : dim; 494 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 495 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 496 ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 497 if (periodicMap) { /* we need to localize coordinates on cells */ 498 ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 499 } else { 500 ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 501 } 502 for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 503 ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 504 ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 505 } 506 if (periodicMap) { 507 ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 508 for (cell = 0, c = 0; c < numCells; ++c) { 509 if (gmsh_elem[c].dim == dim) { 510 PetscInt corner; 511 PetscBool pc = PETSC_FALSE; 512 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 513 pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - shift)); 514 } 515 if (pc) { 516 PetscInt dof = gmsh_elem[c].numNodes*embedDim; 517 ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr); 518 ierr = PetscSectionSetDof(coordSection, cell, dof);CHKERRQ(ierr); 519 ierr = PetscSectionSetFieldDof(coordSection, cell, 0, dof);CHKERRQ(ierr); 520 } 521 cell++; 522 } 523 } 524 } 525 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 526 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 527 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 528 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 529 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 530 ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 531 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 532 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 533 if (periodicMap) { 534 PetscInt off; 535 536 for (cell = 0, c = 0; c < numCells; ++c) { 537 PetscInt pcone[8], corner; 538 if (gmsh_elem[c].dim == dim) { 539 if (PetscUnlikely(PetscBTLookup(periodicC, cell))) { 540 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 541 pcone[corner] = gmsh_elem[c].nodes[corner] - shift; 542 } 543 if (dim == 3) { 544 /* Tetrahedra are inverted */ 545 if (gmsh_elem[c].numNodes == 4) { 546 PetscInt tmp = pcone[0]; 547 pcone[0] = pcone[1]; 548 pcone[1] = tmp; 549 } 550 /* Hexahedra are inverted */ 551 if (gmsh_elem[c].numNodes == 8) { 552 PetscInt tmp = pcone[1]; 553 pcone[1] = pcone[3]; 554 pcone[3] = tmp; 555 } 556 } 557 ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr); 558 for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 559 v = pcone[corner]; 560 for (d = 0; d < embedDim; ++d) { 561 coords[off++] = (PetscReal) coordsIn[v*3+d]; 562 } 563 } 564 } 565 cell++; 566 } 567 } 568 for (v = 0; v < numVertices; ++v) { 569 ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 570 for (d = 0; d < embedDim; ++d) { 571 coords[off+d] = (PetscReal) coordsIn[periodicMapI[v]*3+d]; 572 } 573 } 574 } else { 575 for (v = 0; v < numVertices; ++v) { 576 for (d = 0; d < embedDim; ++d) { 577 coords[v*embedDim+d] = (PetscReal) coordsIn[v*3+d]; 578 } 579 } 580 } 581 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 582 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 583 ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 584 585 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 586 ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 587 ierr = PetscFree(periodicMap);CHKERRQ(ierr); 588 ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 589 ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 590 ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 591 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 592 593 ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 594 PetscFunctionReturn(0); 595 } 596