1 #define PETSCDM_DLL 2 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 4 #undef __FUNCT__ 5 #define __FUNCT__ "DMPlexCreateGmsh" 6 /*@ 7 DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file. 8 9 Collective on comm 10 11 Input Parameters: 12 + comm - The MPI communicator 13 . viewer - The Viewer associated with a Gmsh file 14 - interpolate - Create faces and edges in the mesh 15 16 Output Parameter: 17 . dm - The DM object representing the mesh 18 19 Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 20 21 Level: beginner 22 23 .keywords: mesh,Gmsh 24 .seealso: DMPLEX, DMCreate() 25 @*/ 26 PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 27 { 28 FILE *fd; 29 PetscSection coordSection; 30 Vec coordinates; 31 PetscScalar *coords, *coordsIn = NULL; 32 PetscInt dim = 0, tdim = 0, coordSize, c, v, d, numCorners, cellNum; 33 int numVertices = 0, numCells = 0, trueNumCells = 0, numFacets = 0, cone[8], snum; 34 long fpos = 0; 35 PetscMPIInt num_proc, rank; 36 char line[PETSC_MAX_PATH_LEN]; 37 PetscBool match; 38 PetscErrorCode ierr; 39 40 PetscFunctionBegin; 41 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 42 ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 43 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 44 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 45 if (!rank) { 46 PetscBool match; 47 int fileType, dataSize; 48 49 ierr = PetscViewerASCIIGetPointer(viewer, &fd);CHKERRQ(ierr); 50 /* Read header */ 51 fgets(line, PETSC_MAX_PATH_LEN, fd); 52 ierr = PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 53 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 54 snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2); 55 if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 56 if (dataSize != sizeof(double)) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Data size %d is not valid for a Gmsh file", dataSize); 57 fgets(line, PETSC_MAX_PATH_LEN, fd); 58 ierr = PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 59 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 60 /* Read vertices */ 61 fgets(line, PETSC_MAX_PATH_LEN, fd); 62 ierr = PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 63 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 64 snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1); 65 ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr); 66 for (v = 0; v < numVertices; ++v) { 67 double x, y, z; 68 int i; 69 70 snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4); 71 coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z; 72 if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 73 } 74 fgets(line, PETSC_MAX_PATH_LEN, fd); 75 ierr = PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 76 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 77 /* Read cells */ 78 fgets(line, PETSC_MAX_PATH_LEN, fd); 79 ierr = PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 80 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 81 snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1); 82 } 83 84 if (!rank) { 85 fpos = ftell(fd); 86 /* The Gmsh format disguises facets as elements, so we have to run through all "element" entries 87 to get the correct numCells and decide the topological dimension of the mesh */ 88 trueNumCells = 0; 89 for (c = 0; c < numCells; ++c) { 90 ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone);CHKERRQ(ierr); 91 if (dim > tdim) { 92 tdim = dim; 93 trueNumCells = 0; 94 } 95 trueNumCells++; 96 } 97 } 98 ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 99 numFacets = numCells - trueNumCells; 100 if (!rank) { 101 ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 102 for (c = 0; c < numCells; ++c) { 103 ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone);CHKERRQ(ierr); 104 if (dim == tdim) { 105 ierr = DMPlexSetConeSize(*dm, c-numFacets, numCorners);CHKERRQ(ierr); 106 } 107 if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1); 108 } 109 } 110 ierr = DMSetUp(*dm);CHKERRQ(ierr); 111 if (!rank) { 112 ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 113 PetscInt pcone[8], corner; 114 for (c = 0; c < numCells; ++c) { 115 ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone);CHKERRQ(ierr); 116 if (dim == tdim) { 117 for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + trueNumCells-1; 118 ierr = DMPlexSetCone(*dm, c-numFacets, pcone);CHKERRQ(ierr); 119 } 120 if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1); 121 } 122 fgets(line, PETSC_MAX_PATH_LEN, fd); 123 ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 124 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 125 } 126 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 127 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 128 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 129 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 130 if (interpolate) { 131 DM idm; 132 133 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 134 ierr = DMDestroy(dm);CHKERRQ(ierr); 135 *dm = idm; 136 } 137 /* Read coordinates */ 138 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 139 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 140 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 141 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 142 for (v = numCells; v < numCells+numVertices; ++v) { 143 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 144 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 145 } 146 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 147 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 148 ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 149 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 150 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 151 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 152 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 153 if (!rank) { 154 for (v = 0; v < numVertices; ++v) { 155 for (d = 0; d < dim; ++d) { 156 coords[v*dim+d] = coordsIn[v*3+d]; 157 } 158 } 159 } 160 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 161 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 162 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 163 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 164 PetscFunctionReturn(0); 165 } 166 167 #undef __FUNCT__ 168 #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 169 PetscErrorCode DMPlexCreateGmsh_ReadElement(FILE *fd, PetscInt *dim, PetscInt *cellNum, PetscInt *numCorners, int cone[]) 170 { 171 PetscInt t; 172 int snum, cellType, numTags, tag; 173 174 PetscFunctionBegin; 175 snum = fscanf(fd, "%d %d %d", cellNum, &cellType, &numTags);CHKERRQ(snum != 3); 176 if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);} 177 switch (cellType) { 178 case 1: /* 2-node line */ 179 *dim = 1; 180 *numCorners = 2; 181 snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != *numCorners); 182 break; 183 case 2: /* 3-node triangle */ 184 *dim = 2; 185 *numCorners = 3; 186 snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != *numCorners); 187 break; 188 case 3: /* 4-node quadrangle */ 189 *dim = 2; 190 *numCorners = 4; 191 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners); 192 break; 193 case 4: /* 4-node tetrahedron */ 194 *dim = 3; 195 *numCorners = 4; 196 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners); 197 break; 198 case 5: /* 8-node hexahedron */ 199 *dim = 3; 200 *numCorners = 8; 201 snum = fscanf(fd, "%d %d %d %d %d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3], &cone[4], &cone[5], &cone[6], &cone[7]);CHKERRQ(snum != *numCorners); 202 break; 203 default: 204 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 205 } 206 PetscFunctionReturn(0); 207 } 208