1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2331830f3SMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3331830f3SMatthew G. Knepley 4331830f3SMatthew G. Knepley #undef __FUNCT__ 5331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh" 6331830f3SMatthew G. Knepley /*@ 7331830f3SMatthew G. Knepley DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file. 8331830f3SMatthew G. Knepley 9331830f3SMatthew G. Knepley Collective on comm 10331830f3SMatthew G. Knepley 11331830f3SMatthew G. Knepley Input Parameters: 12331830f3SMatthew G. Knepley + comm - The MPI communicator 13331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 14331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 15331830f3SMatthew G. Knepley 16331830f3SMatthew G. Knepley Output Parameter: 17331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 18331830f3SMatthew G. Knepley 19331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 20331830f3SMatthew G. Knepley 21331830f3SMatthew G. Knepley Level: beginner 22331830f3SMatthew G. Knepley 23331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 24331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 25331830f3SMatthew G. Knepley @*/ 26331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 27331830f3SMatthew G. Knepley { 28331830f3SMatthew G. Knepley FILE *fd; 29331830f3SMatthew G. Knepley PetscSection coordSection; 30331830f3SMatthew G. Knepley Vec coordinates; 31331830f3SMatthew G. Knepley PetscScalar *coords, *coordsIn = NULL; 32331830f3SMatthew G. Knepley PetscInt dim = 0, coordSize, c, v, d; 33331830f3SMatthew G. Knepley int numVertices = 0, numCells = 0, snum; 34331830f3SMatthew G. Knepley long fpos; 35331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 36331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 37331830f3SMatthew G. Knepley PetscBool match; 38331830f3SMatthew G. Knepley PetscErrorCode ierr; 39331830f3SMatthew G. Knepley 40331830f3SMatthew G. Knepley PetscFunctionBegin; 41331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 42331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 43331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 44331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 45331830f3SMatthew G. Knepley if (!rank) { 46331830f3SMatthew G. Knepley PetscBool match; 47331830f3SMatthew G. Knepley int fileType, dataSize; 48331830f3SMatthew G. Knepley 49331830f3SMatthew G. Knepley ierr = PetscViewerASCIIGetPointer(viewer, &fd);CHKERRQ(ierr); 50331830f3SMatthew G. Knepley /* Read header */ 51331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 52331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$MeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 53331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 54331830f3SMatthew G. Knepley snum = fscanf(fd, "2.2 %d %d\n", &fileType, &dataSize);CHKERRQ(snum != 2); 55331830f3SMatthew G. Knepley if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 56331830f3SMatthew 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); 57331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 58331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$EndMeshFormat\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 59331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 60331830f3SMatthew G. Knepley /* Read vertices */ 61331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 62331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$Nodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 63331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 64331830f3SMatthew G. Knepley snum = fscanf(fd, "%d\n", &numVertices);CHKERRQ(snum != 1); 65331830f3SMatthew G. Knepley ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr); 66331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 67331830f3SMatthew G. Knepley double x, y, z; 68331830f3SMatthew G. Knepley int i; 69331830f3SMatthew G. Knepley 70331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %lg %lg %lg\n", &i, &x, &y, &z);CHKERRQ(snum != 4); 71331830f3SMatthew G. Knepley coordsIn[v*3+0] = x; coordsIn[v*3+1] = y; coordsIn[v*3+2] = z; 72331830f3SMatthew G. Knepley if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 73331830f3SMatthew G. Knepley } 74331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 75331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$EndNodes\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 76331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 77331830f3SMatthew G. Knepley /* Read cells */ 78331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 79331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$Elements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 80331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 81331830f3SMatthew G. Knepley snum = fscanf(fd, "%d\n", &numCells);CHKERRQ(snum != 1); 82331830f3SMatthew G. Knepley } 83331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 84331830f3SMatthew G. Knepley if (!rank) { 85331830f3SMatthew G. Knepley fpos = ftell(fd); 86331830f3SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 87331830f3SMatthew G. Knepley PetscInt cone[8], numCorners, t; 88331830f3SMatthew G. Knepley int i, cellType, numTags, tag; 89331830f3SMatthew G. Knepley 90331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3); 91331830f3SMatthew G. Knepley if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);} 92331830f3SMatthew G. Knepley switch (cellType) { 93331830f3SMatthew G. Knepley case 1: /* 2-node line */ 94331830f3SMatthew G. Knepley dim = 1; 95331830f3SMatthew G. Knepley numCorners = 2; 96331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners); 97331830f3SMatthew G. Knepley break; 98331830f3SMatthew G. Knepley case 2: /* 3-node triangle */ 99331830f3SMatthew G. Knepley dim = 2; 100331830f3SMatthew G. Knepley numCorners = 3; 101331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners); 102331830f3SMatthew G. Knepley break; 103331830f3SMatthew G. Knepley case 3: /* 4-node quadrangle */ 104331830f3SMatthew G. Knepley dim = 2; 105331830f3SMatthew G. Knepley numCorners = 4; 106331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 107331830f3SMatthew G. Knepley break; 108331830f3SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 109331830f3SMatthew G. Knepley dim = 3; 110331830f3SMatthew G. Knepley numCorners = 4; 111331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 112331830f3SMatthew G. Knepley break; 113331830f3SMatthew G. Knepley case 5: /* 8-node hexahedron */ 114331830f3SMatthew G. Knepley dim = 3; 115331830f3SMatthew G. Knepley numCorners = 8; 116331830f3SMatthew 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); 117331830f3SMatthew G. Knepley break; 118331830f3SMatthew G. Knepley default: 119331830f3SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 120331830f3SMatthew G. Knepley } 121331830f3SMatthew G. Knepley ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 122331830f3SMatthew G. Knepley if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1); 123331830f3SMatthew G. Knepley } 124331830f3SMatthew G. Knepley } 125331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 126331830f3SMatthew G. Knepley if (!rank) { 127331830f3SMatthew G. Knepley ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 128331830f3SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 129331830f3SMatthew G. Knepley PetscInt cone[8], numCorners, corner, t; 130331830f3SMatthew G. Knepley int i, cellType, numTags, tag; 131331830f3SMatthew G. Knepley 132331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d", &i, &cellType, &numTags);CHKERRQ(snum != 3); 133331830f3SMatthew G. Knepley if (numTags) for (t = 0; t < numTags; ++t) {snum = fscanf(fd, "%d", &tag);CHKERRQ(snum != 1);} 134331830f3SMatthew G. Knepley switch (cellType) { 135331830f3SMatthew G. Knepley case 1: /* 2-node line */ 136331830f3SMatthew G. Knepley dim = 1; 137331830f3SMatthew G. Knepley numCorners = 2; 138331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != numCorners); 139331830f3SMatthew G. Knepley break; 140331830f3SMatthew G. Knepley case 2: /* 3-node triangle */ 141331830f3SMatthew G. Knepley dim = 2; 142331830f3SMatthew G. Knepley numCorners = 3; 143331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != numCorners); 144331830f3SMatthew G. Knepley break; 145331830f3SMatthew G. Knepley case 3: /* 4-node quadrangle */ 146331830f3SMatthew G. Knepley dim = 2; 147331830f3SMatthew G. Knepley numCorners = 4; 148331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 149331830f3SMatthew G. Knepley break; 150331830f3SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 151331830f3SMatthew G. Knepley dim = 3; 152331830f3SMatthew G. Knepley numCorners = 4; 153331830f3SMatthew G. Knepley snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != numCorners); 154331830f3SMatthew G. Knepley break; 155331830f3SMatthew G. Knepley case 5: /* 8-node hexahedron */ 156331830f3SMatthew G. Knepley dim = 3; 157331830f3SMatthew G. Knepley numCorners = 8; 158331830f3SMatthew 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); 159331830f3SMatthew G. Knepley break; 160331830f3SMatthew G. Knepley default: 161331830f3SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 162331830f3SMatthew G. Knepley } 163331830f3SMatthew G. Knepley for (corner = 0; corner < numCorners; ++corner) cone[corner] += numCells-1; 164331830f3SMatthew G. Knepley ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 165331830f3SMatthew G. Knepley if (i != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", i, c+1); 166331830f3SMatthew G. Knepley } 167331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 168331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 169331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 170331830f3SMatthew G. Knepley } 171331830f3SMatthew G. Knepley ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 172331830f3SMatthew G. Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 173331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 174331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 175331830f3SMatthew G. Knepley if (interpolate) { 176331830f3SMatthew G. Knepley DM idm; 177331830f3SMatthew G. Knepley 178331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 179331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 180331830f3SMatthew G. Knepley *dm = idm; 181331830f3SMatthew G. Knepley } 182331830f3SMatthew G. Knepley /* Read coordinates */ 183*979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 184331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 185331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 186331830f3SMatthew G. Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 187331830f3SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) { 188331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 189331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 190331830f3SMatthew G. Knepley } 191331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 192331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 193331830f3SMatthew G. Knepley ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 194331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 195331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 196331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 197331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 198331830f3SMatthew G. Knepley if (!rank) { 199331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 200331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 201331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 202331830f3SMatthew G. Knepley } 203331830f3SMatthew G. Knepley } 204331830f3SMatthew G. Knepley } 205331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 206331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 207331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 208331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 209331830f3SMatthew G. Knepley PetscFunctionReturn(0); 210331830f3SMatthew G. Knepley } 211