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