/*@C DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file + comm - The MPI communicator . filename - Name of the Gmsh file - interpolate - Create faces and edges in the mesh Output Parameter: . dm - The DM object representing the mesh Level: beginner .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() @*/ PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) { PetscViewer viewer, vheader; PetscMPIInt rank; PetscViewerType vtype; char line[PETSC_MAX_PATH_LEN]; int snum; PetscBool match; int fT; PetscInt fileType; float version; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); /* Determine Gmsh file type (ASCII or binary) from file header */ ierr = PetscViewerCreate(comm, &vheader); CHKERRQ(ierr); ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII); CHKERRQ(ierr); ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ); CHKERRQ(ierr); ierr = PetscViewerFileSetName(vheader, filename); CHKERRQ(ierr); if (!rank) { /* Read only the first two lines of the Gmsh file */ ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING); CHKERRQ(ierr); ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match); CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING); CHKERRQ(ierr); snum = sscanf(line, "%f %d", &version, &fT); fileType = (PetscInt) fT; if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); } ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm); CHKERRQ(ierr); /* Create appropriate viewer and build plex */ if (fileType == 0) vtype = PETSCVIEWERASCII; else vtype = PETSCVIEWERBINARY; ierr = PetscViewerCreate(comm, &viewer); CHKERRQ(ierr); ierr = PetscViewerSetType(viewer, vtype); CHKERRQ(ierr); ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ); CHKERRQ(ierr); ierr = PetscViewerFileSetName(viewer, filename); CHKERRQ(ierr); ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm); CHKERRQ(ierr); ierr = PetscViewerDestroy(&viewer); CHKERRQ(ierr); ierr = PetscViewerDestroy(&vheader); CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) { PetscInt ret, i = 0; PetscErrorCode ierr; PetscFunctionBegin; do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR);CHKERRQ(ierr);} while (ret > 0 && buffer[i-1] != '\0' && buffer[i-1] != delim); buffer[i] = '\0'; PetscFunctionReturn(0); }
PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary) { int fdes=0; FILE *file; PetscInt i; PetscErrorCode ierr; PetscFunctionBegin; if (binary) { /* Extract raw file descriptor to read binary block */ ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr); fflush(file); fdes = fileno(file); } if (!binary && dtype == PETSC_INT) { char cbuf[256]; int ibuf, snum; /* Parse hexadecimal ascii integers */ for (i = 0; i < count; i++) { ierr = PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING);CHKERRQ(ierr); snum = sscanf(cbuf, "%x", &ibuf); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); ((PetscInt*)data)[i] = (PetscInt)ibuf; } } else if (binary && dtype == PETSC_INT) { /* Always read 32-bit ints and cast to PetscInt */ int *ibuf; ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr); ierr = PetscBinaryRead(fdes, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr); for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]); ierr = PetscFree(ibuf);CHKERRQ(ierr); } else if (binary && dtype == PETSC_SCALAR) { float *fbuf; /* Always read 32-bit floats and cast to PetscScalar */ ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr); ierr = PetscBinaryRead(fdes, fbuf, count, PETSC_FLOAT);CHKERRQ(ierr); ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr); for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]); ierr = PetscFree(fbuf);CHKERRQ(ierr); } else { ierr = PetscViewerASCIIRead(viewer, data, count, NULL, dtype);CHKERRQ(ierr); } PetscFunctionReturn(0); }
PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) { char buffer[PETSC_MAX_PATH_LEN]; int snum; PetscErrorCode ierr; PetscFunctionBegin; /* Fast-forward to next section and derive its index */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr); snum = sscanf(buffer, "%d", &(s->index)); /* If we can't match an index return -1 to signal end-of-file */ if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} if (s->index == 0) { /* Comment */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } else if (s->index == 2) { /* Dimension */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); snum = sscanf(buffer, "%d", &(s->nd)); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); } else if (s->index == 10 || s->index == 2010) { /* Vertices */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); if (s->zoneID > 0) { PetscInt numCoords = s->last - s->first + 1; ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } else if (s->index == 12 || s->index == 2012) { /* Cells */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); snum = sscanf(buffer, "(%x", &(s->zoneID)); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); if (s->zoneID == 0) { /* Header section */ snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); } else { /* Data section */ snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); if (s->nd == 0) { /* Read cell type definitions for mixed cells */ PetscInt numCells = s->last - s->first + 1; ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr); ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); ierr = PetscFree(s->data);CHKERRQ(ierr); ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } } ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } else if (s->index == 13 || s->index == 2013) { /* Faces */ ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); snum = sscanf(buffer, "(%x", &(s->zoneID)); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); if (s->zoneID == 0) { /* Header section */ snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); } else { /* Data section */ PetscInt f, numEntries, numFaces, numFaceVert; snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); switch (s->nd) { case 0: numEntries = PETSC_DETERMINE; break; case 2: numEntries = 2 + 2; break; /* linear */ case 3: numEntries = 2 + 3; break; /* triangular */ case 4: numEntries = 2 + 4; break; /* quadrilateral */ default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); } numFaces = s->last-s->first + 1; if (numEntries != PETSC_DETERMINE) { /* Allocate space only if we already know the size of the block */ ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); } for (f = 0; f < numFaces; f++) { if (s->nd == 0) { /* Determine the size of the block for "mixed" facets */ ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); if (numEntries == PETSC_DETERMINE) { numEntries = numFaceVert + 2; ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); } else { if (numEntries != numFaceVert + 2) { SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); } } } ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); } s->nd = numEntries - 2; ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); } else { /* Unknown section type */ PetscInt depth = 1; do { /* Match parentheses when parsing unknown sections */ do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);} while (buffer[0] != '(' && buffer[0] != ')'); if (buffer[0] == '(') depth++; if (buffer[0] == ')') depth--; } while (depth > 0); ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); } PetscFunctionReturn(0); }
/*@ DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer Collective on comm Input Parameters: + comm - The MPI communicator . viewer - The Viewer associated with a Gmsh file - interpolate - Create faces and edges in the mesh Output Parameter: . dm - The DM object representing the mesh Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format Level: beginner .keywords: mesh,Gmsh .seealso: DMPLEX, DMCreate() @*/ PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) { PetscViewerType vtype; GmshElement *gmsh_elem; PetscSection coordSection; Vec coordinates; PetscScalar *coords, *coordsIn = NULL; PetscInt dim = 0, coordSize, c, v, d, r, cell; int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum; PetscMPIInt num_proc, rank; char line[PETSC_MAX_PATH_LEN]; PetscBool match, binary, bswap = PETSC_FALSE; PetscErrorCode ierr; PetscFunctionBegin; ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); ierr = DMCreate(comm, dm);CHKERRQ(ierr); ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); if (!rank || binary) { PetscBool match; int fileType, dataSize; float version; /* Read header */ ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); if (binary) { int checkInt; ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); if (checkInt != 1) { ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); if (checkInt == 1) bswap = PETSC_TRUE; else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); } } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); /* OPTIONAL Read physical names */ ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (match) { ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); snum = sscanf(line, "%d", &numRegions); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); for (r = 0; r < numRegions; ++r) { PetscInt rdim, tag; ierr = PetscViewerRead(viewer, &rdim, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); ierr = PetscViewerRead(viewer, &tag, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); } ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); /* Initial read for vertex section */ ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); } /* Read vertices */ ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); snum = sscanf(line, "%d", &numVertices); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); if (binary) { size_t doubleSize, intSize; PetscInt elementSize; char *buffer; PetscScalar *baseptr; ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); elementSize = (intSize + 3*doubleSize); ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); for (v = 0; v < numVertices; ++v) { baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize)); coordsIn[v*3+0] = baseptr[0]; coordsIn[v*3+1] = baseptr[1]; coordsIn[v*3+2] = baseptr[2]; } ierr = PetscFree(buffer);CHKERRQ(ierr); } else { for (v = 0; v < numVertices; ++v) { ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); } } ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); /* Read cells */ ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; snum = sscanf(line, "%d", &numCells); if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); } if (!rank || binary) { /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the file contents multiple times to figure out the true number of cells and facets in the given mesh. To make this more efficient we read the file contents only once and store them in memory, while determining the true number of cells. */ ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); for (trueNumCells=0, c = 0; c < numCells; ++c) { if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} if (gmsh_elem[c].dim == dim) trueNumCells++; } ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); } /* For binary we read on all ranks, but only build the plex on rank 0 */ if (binary && rank) {trueNumCells = 0; numVertices = 0;}; /* Allocate the cell-vertex mesh */ ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); if (!rank) { for (cell = 0, c = 0; c < numCells; ++c) { if (gmsh_elem[c].dim == dim) { ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); cell++; } } } ierr = DMSetUp(*dm);CHKERRQ(ierr); /* Add cell-vertex connections */ if (!rank) { PetscInt pcone[8], corner; for (cell = 0, c = 0; c < numCells; ++c) { if (gmsh_elem[c].dim == dim) { for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1; } if (dim == 3) { /* Tetrahedra are inverted */ if (gmsh_elem[c].numNodes == 4) { PetscInt tmp = pcone[0]; pcone[0] = pcone[1]; pcone[1] = tmp; } /* Hexahedra are inverted */ if (gmsh_elem[c].numNodes == 8) { PetscInt tmp = pcone[1]; pcone[1] = pcone[3]; pcone[3] = tmp; } } ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); cell++; } } } ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); ierr = DMPlexStratify(*dm);CHKERRQ(ierr); if (interpolate) { DM idm = NULL; ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); ierr = DMDestroy(dm);CHKERRQ(ierr); *dm = idm; } if (!rank) { /* Apply boundary IDs by finding the relevant facets with vertex joins */ PetscInt pcone[8], corner, vStart, vEnd; ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); for (c = 0; c < numCells; ++c) { if (gmsh_elem[c].dim == dim-1) { PetscInt joinSize; const PetscInt *join; for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1; } ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); } } /* Create cell sets */ for (cell = 0, c = 0; c < numCells; ++c) { if (gmsh_elem[c].dim == dim) { if (gmsh_elem[c].numTags > 0) { ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); cell++; } } } } /* Read coordinates */ ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); } ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); if (!rank) { for (v = 0; v < numVertices; ++v) { for (d = 0; d < dim; ++d) { coords[v*dim+d] = coordsIn[v*3+d]; } } } ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); ierr = PetscFree(coordsIn);CHKERRQ(ierr); ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); ierr = VecDestroy(&coordinates);CHKERRQ(ierr); /* Clean up intermediate storage */ if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); PetscFunctionReturn(0); }
PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems) { PetscInt c, p; GmshElement *elements; int i, cellType, dim, numNodes, numElem, numTags; int ibuf[16]; PetscErrorCode ierr; PetscFunctionBegin; ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); for (c = 0; c < numCells;) { ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr); if (binary) { cellType = ibuf[0]; numElem = ibuf[1]; numTags = ibuf[2]; } else { elements[c].id = ibuf[0]; cellType = ibuf[1]; numTags = ibuf[2]; numElem = 1; } switch (cellType) { case 1: /* 2-node line */ dim = 1; numNodes = 2; break; case 2: /* 3-node triangle */ dim = 2; numNodes = 3; break; case 3: /* 4-node quadrangle */ dim = 2; numNodes = 4; break; case 4: /* 4-node tetrahedron */ dim = 3; numNodes = 4; break; case 5: /* 8-node hexahedron */ dim = 3; numNodes = 8; break; case 15: /* 1-node vertex */ dim = 0; numNodes = 1; break; default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); } if (binary) { const PetscInt nint = numNodes + numTags + 1; for (i = 0; i < numElem; ++i, ++c) { /* Loop over inner binary element block */ elements[c].dim = dim; elements[c].numNodes = numNodes; elements[c].numTags = numTags; ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr); elements[c].id = ibuf[0]; for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; } } else { elements[c].dim = dim; elements[c].numNodes = numNodes; elements[c].numTags = numTags; ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); c++; } } *gmsh_elems = elements; PetscFunctionReturn(0); }