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, coordSize, c, v, d; 33 int numVertices = 0, numCells = 0, 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 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 84 if (!rank) { 85 fpos = ftell(fd); 86 for (c = 0; c < numCells; ++c) { 87 PetscInt numCorners, t; 88 int cone[8], i, cellType, numTags, tag; 89 90 snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3); 91 if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);} 92 switch (cellType) { 93 case 1: /* 2-node line */ 94 dim = 1; 95 numCorners = 2; 96 snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners); 97 break; 98 case 2: /* 3-node triangle */ 99 dim = 2; 100 numCorners = 3; 101 snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners); 102 break; 103 case 3: /* 4-node quadrangle */ 104 dim = 2; 105 numCorners = 4; 106 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 107 break; 108 case 4: /* 4-node tetrahedron */ 109 dim = 3; 110 numCorners = 4; 111 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 112 break; 113 case 5: /* 8-node hexahedron */ 114 dim = 3; 115 numCorners = 8; 116 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); 117 break; 118 default: 119 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 120 } 121 ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 122 if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1); 123 } 124 } 125 ierr = DMSetUp(*dm);CHKERRQ(ierr); 126 if (!rank) { 127 ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 128 for (c = 0; c < numCells; ++c) { 129 PetscInt pcone[8], numCorners, corner, t; 130 int cone[8], i, cellType, numTags, tag; 131 132 snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3); 133 if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);} 134 switch (cellType) { 135 case 1: /* 2-node line */ 136 dim = 1; 137 numCorners = 2; 138 snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners); 139 break; 140 case 2: /* 3-node triangle */ 141 dim = 2; 142 numCorners = 3; 143 snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners); 144 break; 145 case 3: /* 4-node quadrangle */ 146 dim = 2; 147 numCorners = 4; 148 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 149 break; 150 case 4: /* 4-node tetrahedron */ 151 dim = 3; 152 numCorners = 4; 153 snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 154 ierr = DMPlexInvertCell(dim, numCorners, cone);CHKERRQ(ierr); 155 break; 156 case 5: /* 8-node hexahedron */ 157 dim = 3; 158 numCorners = 8; 159 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); 160 ierr = DMPlexInvertCell(dim, numCorners, cone);CHKERRQ(ierr); 161 break; 162 default: 163 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 164 } 165 for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + numCells-1; 166 ierr = DMPlexSetCone(*dm, c, pcone);CHKERRQ(ierr); 167 if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1); 168 } 169 fgets(line, PETSC_MAX_PATH_LEN, fd); 170 ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 171 if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 172 } 173 ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 174 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 175 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 176 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 177 if (interpolate) { 178 DM idm; 179 180 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 181 ierr = DMDestroy(dm);CHKERRQ(ierr); 182 *dm = idm; 183 } 184 /* Read coordinates */ 185 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 186 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 187 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 188 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 189 for (v = numCells; v < numCells+numVertices; ++v) { 190 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 191 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 192 } 193 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 194 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 195 ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 196 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 197 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 198 ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 199 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 200 if (!rank) { 201 for (v = 0; v < numVertices; ++v) { 202 for (d = 0; d < dim; ++d) { 203 coords[v*dim+d] = coordsIn[v*3+d]; 204 } 205 } 206 } 207 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 208 ierr = PetscFree(coordsIn);CHKERRQ(ierr); 209 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 210 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 211 PetscFunctionReturn(0); 212 } 213