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; 32*16ad7e67SMichael Lange PetscInt dim = 0, tdim = 0, coordSize, c, v, d, cellNum, numCorners, numTags; 33*16ad7e67SMichael Lange int numVertices = 0, numCells = 0, trueNumCells = 0, numFacets = 0, cone[8], tags[2], snum; 34261b4668SMatthew G. Knepley long fpos = 0; 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 } 83db397164SMichael Lange 84331830f3SMatthew G. Knepley if (!rank) { 85331830f3SMatthew G. Knepley fpos = ftell(fd); 86db397164SMichael Lange /* The Gmsh format disguises facets as elements, so we have to run through all "element" entries 87db397164SMichael Lange to get the correct numCells and decide the topological dimension of the mesh */ 88db397164SMichael Lange trueNumCells = 0; 89db397164SMichael Lange for (c = 0; c < numCells; ++c) { 90*16ad7e67SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, &numTags, tags);CHKERRQ(ierr); 91db397164SMichael Lange if (dim > tdim) { 92db397164SMichael Lange tdim = dim; 93db397164SMichael Lange trueNumCells = 0; 94db397164SMichael Lange } 95db397164SMichael Lange trueNumCells++; 96db397164SMichael Lange } 97db397164SMichael Lange } 98db397164SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 99db397164SMichael Lange numFacets = numCells - trueNumCells; 100db397164SMichael Lange if (!rank) { 101db397164SMichael Lange ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 102331830f3SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 103*16ad7e67SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, &numTags, tags);CHKERRQ(ierr); 104db397164SMichael Lange if (dim == tdim) { 105db397164SMichael Lange ierr = DMPlexSetConeSize(*dm, c-numFacets, numCorners);CHKERRQ(ierr); 106db397164SMichael Lange } 10730412aabSMichael Lange if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1); 108331830f3SMatthew G. Knepley } 109331830f3SMatthew G. Knepley } 110331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 111331830f3SMatthew G. Knepley if (!rank) { 112331830f3SMatthew G. Knepley ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 11330412aabSMichael Lange PetscInt pcone[8], corner; 114331830f3SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 115*16ad7e67SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, &numTags, tags);CHKERRQ(ierr); 116db397164SMichael Lange if (dim == tdim) { 117db397164SMichael Lange for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + trueNumCells-1; 118db397164SMichael Lange ierr = DMPlexSetCone(*dm, c-numFacets, pcone);CHKERRQ(ierr); 119db397164SMichael Lange } 12030412aabSMichael Lange if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1); 121331830f3SMatthew G. Knepley } 122331830f3SMatthew G. Knepley fgets(line, PETSC_MAX_PATH_LEN, fd); 123331830f3SMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 124331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 125331830f3SMatthew G. Knepley } 1266228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 127331830f3SMatthew G. Knepley ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 128331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 129331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 130331830f3SMatthew G. Knepley if (interpolate) { 131331830f3SMatthew G. Knepley DM idm; 132331830f3SMatthew G. Knepley 133331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 134331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 135331830f3SMatthew G. Knepley *dm = idm; 136331830f3SMatthew G. Knepley } 137*16ad7e67SMichael Lange 138*16ad7e67SMichael Lange if (!rank) { 139*16ad7e67SMichael Lange ierr = fseek(fd, fpos, SEEK_SET);CHKERRQ(ierr); 140*16ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 141*16ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 142*16ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 143*16ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 144*16ad7e67SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(fd, &dim, &cellNum, &numCorners, cone, &numTags, tags);CHKERRQ(ierr); 145*16ad7e67SMichael Lange if (dim == tdim-1) { 146*16ad7e67SMichael Lange PetscInt joinSize; 147*16ad7e67SMichael Lange const PetscInt *join; 148*16ad7e67SMichael Lange for (corner = 0; corner < numCorners; ++corner) pcone[corner] = cone[corner] + vStart - 1; 149*16ad7e67SMichael Lange ierr = DMPlexGetFullJoin(*dm, numCorners, pcone, &joinSize, &join);CHKERRQ(ierr); 150*16ad7e67SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", cellNum); 151*16ad7e67SMichael Lange ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], tags[0]);CHKERRQ(ierr); 152*16ad7e67SMichael Lange ierr = DMPlexRestoreJoin(*dm, numCorners, pcone, &joinSize, &join);CHKERRQ(ierr); 153*16ad7e67SMichael Lange } 154*16ad7e67SMichael Lange if (cellNum != c+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell number %d should be %d", cellNum, c+1); 155*16ad7e67SMichael Lange } 156*16ad7e67SMichael Lange fgets(line, PETSC_MAX_PATH_LEN, fd); 157*16ad7e67SMichael Lange ierr = PetscStrncmp(line, "$EndElements\n", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 158*16ad7e67SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 159*16ad7e67SMichael Lange } 160*16ad7e67SMichael Lange 161331830f3SMatthew G. Knepley /* Read coordinates */ 162979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 163331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 164331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 165331830f3SMatthew G. Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 166331830f3SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) { 167331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 168331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 169331830f3SMatthew G. Knepley } 170331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 171331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 172331830f3SMatthew G. Knepley ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 173331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 174331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 175331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 176331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 177331830f3SMatthew G. Knepley if (!rank) { 178331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 179331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 180331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 181331830f3SMatthew G. Knepley } 182331830f3SMatthew G. Knepley } 183331830f3SMatthew G. Knepley } 184331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 185331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 186331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 187331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 188331830f3SMatthew G. Knepley PetscFunctionReturn(0); 189331830f3SMatthew G. Knepley } 190db397164SMichael Lange 191db397164SMichael Lange #undef __FUNCT__ 19230412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 193*16ad7e67SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(FILE *fd, PetscInt *dim, PetscInt *cellNum, PetscInt *numCorners, int cone[], PetscInt *numTags, int tags[]) 194db397164SMichael Lange { 19530412aabSMichael Lange PetscInt t; 196*16ad7e67SMichael Lange int snum, cellType; 197db397164SMichael Lange 19830412aabSMichael Lange PetscFunctionBegin; 199*16ad7e67SMichael Lange snum = fscanf(fd, "%d %d %d", cellNum, &cellType, numTags);CHKERRQ(snum != 3); 200*16ad7e67SMichael Lange if (numTags) for (t = 0; t < *numTags; ++t) {snum = fscanf(fd, "%d", &tags[t]);CHKERRQ(snum != 1);} 20130412aabSMichael Lange switch (cellType) { 20230412aabSMichael Lange case 1: /* 2-node line */ 20330412aabSMichael Lange *dim = 1; 20430412aabSMichael Lange *numCorners = 2; 20530412aabSMichael Lange snum = fscanf(fd, "%d %d\n", &cone[0], &cone[1]);CHKERRQ(snum != *numCorners); 20630412aabSMichael Lange break; 20730412aabSMichael Lange case 2: /* 3-node triangle */ 20830412aabSMichael Lange *dim = 2; 20930412aabSMichael Lange *numCorners = 3; 21030412aabSMichael Lange snum = fscanf(fd, "%d %d %d\n", &cone[0], &cone[1], &cone[2]);CHKERRQ(snum != *numCorners); 21130412aabSMichael Lange break; 21230412aabSMichael Lange case 3: /* 4-node quadrangle */ 21330412aabSMichael Lange *dim = 2; 21430412aabSMichael Lange *numCorners = 4; 21530412aabSMichael Lange snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners); 21630412aabSMichael Lange break; 21730412aabSMichael Lange case 4: /* 4-node tetrahedron */ 21830412aabSMichael Lange *dim = 3; 21930412aabSMichael Lange *numCorners = 4; 22030412aabSMichael Lange snum = fscanf(fd, "%d %d %d %d\n", &cone[0], &cone[1], &cone[2], &cone[3]);CHKERRQ(snum != *numCorners); 22130412aabSMichael Lange break; 22230412aabSMichael Lange case 5: /* 8-node hexahedron */ 22330412aabSMichael Lange *dim = 3; 22430412aabSMichael Lange *numCorners = 8; 22530412aabSMichael Lange 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); 22630412aabSMichael Lange break; 22730412aabSMichael Lange default: 22830412aabSMichael Lange SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 22930412aabSMichael Lange } 23030412aabSMichael Lange PetscFunctionReturn(0); 231db397164SMichael Lange } 232