1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3331830f3SMatthew G. Knepley 47d282ae0SMichael Lange /*@C 57d282ae0SMichael Lange DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 67d282ae0SMichael Lange 77d282ae0SMichael Lange + comm - The MPI communicator 87d282ae0SMichael Lange . filename - Name of the Gmsh file 97d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh 107d282ae0SMichael Lange 117d282ae0SMichael Lange Output Parameter: 127d282ae0SMichael Lange . dm - The DM object representing the mesh 137d282ae0SMichael Lange 147d282ae0SMichael Lange Level: beginner 157d282ae0SMichael Lange 167d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 177d282ae0SMichael Lange @*/ 187d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 197d282ae0SMichael Lange { 207d282ae0SMichael Lange PetscViewer viewer, vheader; 21abc86ac4SMichael Lange PetscMPIInt rank; 227d282ae0SMichael Lange PetscViewerType vtype; 237d282ae0SMichael Lange char line[PETSC_MAX_PATH_LEN]; 247d282ae0SMichael Lange int snum; 257d282ae0SMichael Lange PetscBool match; 26c1c22fd2SMatthew G. Knepley int fT; 27c1c22fd2SMatthew G. Knepley PetscInt fileType; 28f6ac7a6aSMichael Lange float version; 297d282ae0SMichael Lange PetscErrorCode ierr; 307d282ae0SMichael Lange 317d282ae0SMichael Lange PetscFunctionBegin; 32abc86ac4SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 337d282ae0SMichael Lange /* Determine Gmsh file type (ASCII or binary) from file header */ 347d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr); 357d282ae0SMichael Lange ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr); 367d282ae0SMichael Lange ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr); 377d282ae0SMichael Lange ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr); 38abc86ac4SMichael Lange if (!rank) { 397d282ae0SMichael Lange /* Read only the first two lines of the Gmsh file */ 40060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 417d282ae0SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 427d282ae0SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 43060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 44c1c22fd2SMatthew G. Knepley snum = sscanf(line, "%f %d", &version, &fT); 45c1c22fd2SMatthew G. Knepley fileType = (PetscInt) fT; 46f6ac7a6aSMichael Lange if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 47f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 48abc86ac4SMichael Lange } 49abc86ac4SMichael Lange ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 507d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 517d282ae0SMichael Lange if (fileType == 0) vtype = PETSCVIEWERASCII; 527d282ae0SMichael Lange else vtype = PETSCVIEWERBINARY; 537d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 547d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 557d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 567d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 577d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 587d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 607d282ae0SMichael Lange PetscFunctionReturn(0); 617d282ae0SMichael Lange } 627d282ae0SMichael Lange 63*df032b82SMatthew G. Knepley static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems) 64*df032b82SMatthew G. Knepley { 65*df032b82SMatthew G. Knepley PetscInt c, p; 66*df032b82SMatthew G. Knepley GmshElement *elements; 67*df032b82SMatthew G. Knepley int i, cellType, dim, numNodes, numElem, numTags; 68*df032b82SMatthew G. Knepley int ibuf[16]; 69*df032b82SMatthew G. Knepley PetscErrorCode ierr; 70*df032b82SMatthew G. Knepley 71*df032b82SMatthew G. Knepley PetscFunctionBegin; 72*df032b82SMatthew G. Knepley ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 73*df032b82SMatthew G. Knepley for (c = 0; c < numCells;) { 74*df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 75*df032b82SMatthew G. Knepley if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr); 76*df032b82SMatthew G. Knepley if (binary) { 77*df032b82SMatthew G. Knepley cellType = ibuf[0]; 78*df032b82SMatthew G. Knepley numElem = ibuf[1]; 79*df032b82SMatthew G. Knepley numTags = ibuf[2]; 80*df032b82SMatthew G. Knepley } else { 81*df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 82*df032b82SMatthew G. Knepley cellType = ibuf[1]; 83*df032b82SMatthew G. Knepley numTags = ibuf[2]; 84*df032b82SMatthew G. Knepley numElem = 1; 85*df032b82SMatthew G. Knepley } 86*df032b82SMatthew G. Knepley switch (cellType) { 87*df032b82SMatthew G. Knepley case 1: /* 2-node line */ 88*df032b82SMatthew G. Knepley dim = 1; 89*df032b82SMatthew G. Knepley numNodes = 2; 90*df032b82SMatthew G. Knepley break; 91*df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 92*df032b82SMatthew G. Knepley dim = 2; 93*df032b82SMatthew G. Knepley numNodes = 3; 94*df032b82SMatthew G. Knepley break; 95*df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 96*df032b82SMatthew G. Knepley dim = 2; 97*df032b82SMatthew G. Knepley numNodes = 4; 98*df032b82SMatthew G. Knepley break; 99*df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 100*df032b82SMatthew G. Knepley dim = 3; 101*df032b82SMatthew G. Knepley numNodes = 4; 102*df032b82SMatthew G. Knepley break; 103*df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 104*df032b82SMatthew G. Knepley dim = 3; 105*df032b82SMatthew G. Knepley numNodes = 8; 106*df032b82SMatthew G. Knepley break; 107*df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 108*df032b82SMatthew G. Knepley dim = 0; 109*df032b82SMatthew G. Knepley numNodes = 1; 110*df032b82SMatthew G. Knepley break; 111*df032b82SMatthew G. Knepley default: 112*df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 113*df032b82SMatthew G. Knepley } 114*df032b82SMatthew G. Knepley if (binary) { 115*df032b82SMatthew G. Knepley const PetscInt nint = numNodes + numTags + 1; 116*df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 117*df032b82SMatthew G. Knepley /* Loop over inner binary element block */ 118*df032b82SMatthew G. Knepley elements[c].dim = dim; 119*df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 120*df032b82SMatthew G. Knepley elements[c].numTags = numTags; 121*df032b82SMatthew G. Knepley 122*df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 123*df032b82SMatthew G. Knepley if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr); 124*df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 125*df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 126*df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 127*df032b82SMatthew G. Knepley } 128*df032b82SMatthew G. Knepley } else { 129*df032b82SMatthew G. Knepley elements[c].dim = dim; 130*df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 131*df032b82SMatthew G. Knepley elements[c].numTags = numTags; 132*df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 133*df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 134*df032b82SMatthew G. Knepley c++; 135*df032b82SMatthew G. Knepley } 136*df032b82SMatthew G. Knepley } 137*df032b82SMatthew G. Knepley *gmsh_elems = elements; 138*df032b82SMatthew G. Knepley PetscFunctionReturn(0); 139*df032b82SMatthew G. Knepley } 1407d282ae0SMichael Lange 141331830f3SMatthew G. Knepley /*@ 1427d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 143331830f3SMatthew G. Knepley 144331830f3SMatthew G. Knepley Collective on comm 145331830f3SMatthew G. Knepley 146331830f3SMatthew G. Knepley Input Parameters: 147331830f3SMatthew G. Knepley + comm - The MPI communicator 148331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 149331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 150331830f3SMatthew G. Knepley 151331830f3SMatthew G. Knepley Output Parameter: 152331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 153331830f3SMatthew G. Knepley 154331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 1553d51f2daSMichael Lange and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 156331830f3SMatthew G. Knepley 157331830f3SMatthew G. Knepley Level: beginner 158331830f3SMatthew G. Knepley 159331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 160331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 161331830f3SMatthew G. Knepley @*/ 162331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 163331830f3SMatthew G. Knepley { 16404d1ad83SMichael Lange PetscViewerType vtype; 165a4bb7517SMichael Lange GmshElement *gmsh_elem; 166331830f3SMatthew G. Knepley PetscSection coordSection; 167331830f3SMatthew G. Knepley Vec coordinates; 168331830f3SMatthew G. Knepley PetscScalar *coords, *coordsIn = NULL; 169dc0ede3bSMatthew G. Knepley PetscInt dim = 0, coordSize, c, v, d, r, cell; 170dc0ede3bSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum; 171331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 172331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 17304d1ad83SMichael Lange PetscBool match, binary, bswap = PETSC_FALSE; 174331830f3SMatthew G. Knepley PetscErrorCode ierr; 175331830f3SMatthew G. Knepley 176331830f3SMatthew G. Knepley PetscFunctionBegin; 177331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 178331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 179331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 180331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1813b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 18204d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 18304d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 184abc86ac4SMichael Lange if (!rank || binary) { 185331830f3SMatthew G. Knepley PetscBool match; 186331830f3SMatthew G. Knepley int fileType, dataSize; 187f6ac7a6aSMichael Lange float version; 188331830f3SMatthew G. Knepley 189331830f3SMatthew G. Knepley /* Read header */ 190060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 19104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 192331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 193060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 194f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 195f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 196f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 197331830f3SMatthew 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); 19804d1ad83SMichael Lange if (binary) { 199b9eae255SMichael Lange int checkInt; 200060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 20104d1ad83SMichael Lange if (checkInt != 1) { 202b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 20304d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 20404d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 20504d1ad83SMichael Lange } 2060877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 207060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 20804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 209331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 210dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 211060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 212dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 213dc0ede3bSMatthew G. Knepley if (match) { 214dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 215dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 216dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 217a8667449SSander Arens for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);} 218dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 219dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 220dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 221dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 222dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 223dc0ede3bSMatthew G. Knepley } 224dc0ede3bSMatthew G. Knepley /* Read vertices */ 22504d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 226331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 227060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 2280877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 2290877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 230854ce69bSBarry Smith ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); 23113aecfe5SMichael Lange if (binary) { 23213aecfe5SMichael Lange size_t doubleSize, intSize; 23313aecfe5SMichael Lange PetscInt elementSize; 23413aecfe5SMichael Lange char *buffer; 23513aecfe5SMichael Lange PetscScalar *baseptr; 23613aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); 23713aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); 23813aecfe5SMichael Lange elementSize = (intSize + 3*doubleSize); 23913aecfe5SMichael Lange ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); 24013aecfe5SMichael Lange ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); 24113aecfe5SMichael Lange if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); 242331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 24313aecfe5SMichael Lange baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize)); 24413aecfe5SMichael Lange coordsIn[v*3+0] = baseptr[0]; 24513aecfe5SMichael Lange coordsIn[v*3+1] = baseptr[1]; 24613aecfe5SMichael Lange coordsIn[v*3+2] = baseptr[2]; 24713aecfe5SMichael Lange } 24813aecfe5SMichael Lange ierr = PetscFree(buffer);CHKERRQ(ierr); 24913aecfe5SMichael Lange } else { 25013aecfe5SMichael Lange for (v = 0; v < numVertices; ++v) { 251060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 252060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 253b9eae255SMichael Lange if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 254331830f3SMatthew G. Knepley } 25513aecfe5SMichael Lange } 256060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 25704d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 258331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 259331830f3SMatthew G. Knepley /* Read cells */ 260060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 26104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 262331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 263060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 2640877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 2650877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 266331830f3SMatthew G. Knepley } 267db397164SMichael Lange 268abc86ac4SMichael Lange if (!rank || binary) { 269a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 270a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 271a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 272a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 2733d51f2daSMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); 274a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 275a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 276a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 277db397164SMichael Lange } 278060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 2791b82c3ebSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 280331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 281331830f3SMatthew G. Knepley } 282abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 283abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 284a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 285331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 286331830f3SMatthew G. Knepley if (!rank) { 287a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 288a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 289a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 290a4bb7517SMichael Lange cell++; 291331830f3SMatthew G. Knepley } 292331830f3SMatthew G. Knepley } 293331830f3SMatthew G. Knepley } 294331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 295a4bb7517SMichael Lange /* Add cell-vertex connections */ 296331830f3SMatthew G. Knepley if (!rank) { 297331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 298a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 299a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 300a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 301a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1; 302331830f3SMatthew G. Knepley } 30397ed6be6Ssarens if (dim == 3) { 30497ed6be6Ssarens /* Tetrahedra are inverted */ 30597ed6be6Ssarens if (gmsh_elem[c].numNodes == 4) { 30697ed6be6Ssarens PetscInt tmp = pcone[0]; 30797ed6be6Ssarens pcone[0] = pcone[1]; 30897ed6be6Ssarens pcone[1] = tmp; 30997ed6be6Ssarens } 31097ed6be6Ssarens /* Hexahedra are inverted */ 31197ed6be6Ssarens if (gmsh_elem[c].numNodes == 8) { 31297ed6be6Ssarens PetscInt tmp = pcone[1]; 31397ed6be6Ssarens pcone[1] = pcone[3]; 31497ed6be6Ssarens pcone[3] = tmp; 31597ed6be6Ssarens } 31697ed6be6Ssarens } 317a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 318a4bb7517SMichael Lange cell++; 319331830f3SMatthew G. Knepley } 320a4bb7517SMichael Lange } 321331830f3SMatthew G. Knepley } 3226228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 323c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 324331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 325331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 326331830f3SMatthew G. Knepley if (interpolate) { 327e56d480eSMatthew G. Knepley DM idm = NULL; 328331830f3SMatthew G. Knepley 329331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 330331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 331331830f3SMatthew G. Knepley *dm = idm; 332331830f3SMatthew G. Knepley } 33316ad7e67SMichael Lange 33416ad7e67SMichael Lange if (!rank) { 33516ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 33616ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 337d1a54cefSMatthew G. Knepley 33816ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 33916ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 340a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 34116ad7e67SMichael Lange PetscInt joinSize; 34216ad7e67SMichael Lange const PetscInt *join; 343a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 344a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1; 34516ad7e67SMichael Lange } 346a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 347a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 348c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 349a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 35016ad7e67SMichael Lange } 351a4bb7517SMichael Lange } 3526e1dd89aSLawrence Mitchell 3536e1dd89aSLawrence Mitchell /* Create cell sets */ 3546e1dd89aSLawrence Mitchell for (cell = 0, c = 0; c < numCells; ++c) { 3556e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 3566e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 3576e1dd89aSLawrence Mitchell ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 3586e1dd89aSLawrence Mitchell cell++; 3596e1dd89aSLawrence Mitchell } 3606e1dd89aSLawrence Mitchell } 3616e1dd89aSLawrence Mitchell } 3620c070f12SSander Arens 3630c070f12SSander Arens /* Create vertex sets */ 3640c070f12SSander Arens for (c = 0; c < numCells; ++c) { 3650c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 3660c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 3670c070f12SSander Arens ierr = DMSetLabelValue(*dm, "Vertex Sets", gmsh_elem[c].nodes[0] + vStart - 1, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 3680c070f12SSander Arens } 3690c070f12SSander Arens } 3700c070f12SSander Arens } 37116ad7e67SMichael Lange } 37216ad7e67SMichael Lange 373331830f3SMatthew G. Knepley /* Read coordinates */ 374979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 375331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 376331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 37775b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 37875b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 379331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 380331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 381331830f3SMatthew G. Knepley } 382331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 383331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 3848b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 385331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 386331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 3877635fff5Ssarens ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 388331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 389331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 390331830f3SMatthew G. Knepley if (!rank) { 391331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 392331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 393331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 394331830f3SMatthew G. Knepley } 395331830f3SMatthew G. Knepley } 396331830f3SMatthew G. Knepley } 397331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 398331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 399331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 400331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 401a4bb7517SMichael Lange /* Clean up intermediate storage */ 40229403c87SMichael Lange if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 4033b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 404331830f3SMatthew G. Knepley PetscFunctionReturn(0); 405331830f3SMatthew G. Knepley } 406