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