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__ 57d282ae0SMichael Lange #define __FUNCT__ "DMPlexCreateGmshFromFile" 67d282ae0SMichael Lange /*@C 77d282ae0SMichael Lange DMPlexCreateGmshFromFile - Create a DMPlex mesh from a Gmsh file 87d282ae0SMichael Lange 97d282ae0SMichael Lange + comm - The MPI communicator 107d282ae0SMichael Lange . filename - Name of the Gmsh file 117d282ae0SMichael Lange - interpolate - Create faces and edges in the mesh 127d282ae0SMichael Lange 137d282ae0SMichael Lange Output Parameter: 147d282ae0SMichael Lange . dm - The DM object representing the mesh 157d282ae0SMichael Lange 167d282ae0SMichael Lange Level: beginner 177d282ae0SMichael Lange 187d282ae0SMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 197d282ae0SMichael Lange @*/ 207d282ae0SMichael Lange PetscErrorCode DMPlexCreateGmshFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 217d282ae0SMichael Lange { 227d282ae0SMichael Lange PetscViewer viewer, vheader; 23abc86ac4SMichael Lange PetscMPIInt rank; 247d282ae0SMichael Lange PetscViewerType vtype; 257d282ae0SMichael Lange char line[PETSC_MAX_PATH_LEN]; 267d282ae0SMichael Lange int snum; 277d282ae0SMichael Lange PetscBool match; 28c1c22fd2SMatthew G. Knepley int fT; 29c1c22fd2SMatthew G. Knepley PetscInt fileType; 30f6ac7a6aSMichael Lange float version; 317d282ae0SMichael Lange PetscErrorCode ierr; 327d282ae0SMichael Lange 337d282ae0SMichael Lange PetscFunctionBegin; 34abc86ac4SMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 357d282ae0SMichael Lange /* Determine Gmsh file type (ASCII or binary) from file header */ 367d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &vheader);CHKERRQ(ierr); 377d282ae0SMichael Lange ierr = PetscViewerSetType(vheader, PETSCVIEWERASCII);CHKERRQ(ierr); 387d282ae0SMichael Lange ierr = PetscViewerFileSetMode(vheader, FILE_MODE_READ);CHKERRQ(ierr); 397d282ae0SMichael Lange ierr = PetscViewerFileSetName(vheader, filename);CHKERRQ(ierr); 40abc86ac4SMichael Lange if (!rank) { 41*f8e4bde8SMatthew G. Knepley PetscInt n = 1, m = 2; 427d282ae0SMichael Lange /* Read only the first two lines of the Gmsh file */ 43*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, &n, PETSC_STRING);CHKERRQ(ierr); 447d282ae0SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 457d282ae0SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 46*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, &m, PETSC_STRING);CHKERRQ(ierr); 47c1c22fd2SMatthew G. Knepley snum = sscanf(line, "%f %d", &version, &fT); 48c1c22fd2SMatthew G. Knepley fileType = (PetscInt) fT; 49f6ac7a6aSMichael Lange if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 50f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 51abc86ac4SMichael Lange } 52abc86ac4SMichael Lange ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 537d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 547d282ae0SMichael Lange if (fileType == 0) vtype = PETSCVIEWERASCII; 557d282ae0SMichael Lange else vtype = PETSCVIEWERBINARY; 567d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 577d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 587d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 607d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 617d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 627d282ae0SMichael Lange ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 637d282ae0SMichael Lange PetscFunctionReturn(0); 647d282ae0SMichael Lange } 657d282ae0SMichael Lange 667d282ae0SMichael Lange 677d282ae0SMichael Lange #undef __FUNCT__ 68331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh" 69331830f3SMatthew G. Knepley /*@ 707d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 71331830f3SMatthew G. Knepley 72331830f3SMatthew G. Knepley Collective on comm 73331830f3SMatthew G. Knepley 74331830f3SMatthew G. Knepley Input Parameters: 75331830f3SMatthew G. Knepley + comm - The MPI communicator 76331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 77331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 78331830f3SMatthew G. Knepley 79331830f3SMatthew G. Knepley Output Parameter: 80331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 81331830f3SMatthew G. Knepley 82331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 83331830f3SMatthew G. Knepley 84331830f3SMatthew G. Knepley Level: beginner 85331830f3SMatthew G. Knepley 86331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 87331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 88331830f3SMatthew G. Knepley @*/ 89331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 90331830f3SMatthew G. Knepley { 9104d1ad83SMichael Lange PetscViewerType vtype; 92a4bb7517SMichael Lange GmshElement *gmsh_elem; 93331830f3SMatthew G. Knepley PetscSection coordSection; 94331830f3SMatthew G. Knepley Vec coordinates; 95331830f3SMatthew G. Knepley PetscScalar *coords, *coordsIn = NULL; 96a4bb7517SMichael Lange PetscInt dim = 0, coordSize, c, v, d, cell; 97a4bb7517SMichael Lange int numVertices = 0, numCells = 0, trueNumCells = 0, snum; 98331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 99331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 10004d1ad83SMichael Lange PetscBool match, binary, bswap = PETSC_FALSE; 101331830f3SMatthew G. Knepley PetscErrorCode ierr; 102331830f3SMatthew G. Knepley 103331830f3SMatthew G. Knepley PetscFunctionBegin; 104331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 105331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 106331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 107331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 10804d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 10904d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 110abc86ac4SMichael Lange if (!rank || binary) { 111331830f3SMatthew G. Knepley PetscBool match; 112331830f3SMatthew G. Knepley int fileType, dataSize; 113f6ac7a6aSMichael Lange float version; 114*f8e4bde8SMatthew G. Knepley PetscInt n = 1, m = 3; 115331830f3SMatthew G. Knepley 116331830f3SMatthew G. Knepley /* Read header */ 117*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr); 11804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 119331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 120*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &m, PETSC_STRING);CHKERRQ(ierr); 121f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 122f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 123f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 124331830f3SMatthew 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); 12504d1ad83SMichael Lange if (binary) { 126b9eae255SMichael Lange int checkInt; 127*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, &n, PETSC_ENUM);CHKERRQ(ierr); 12804d1ad83SMichael Lange if (checkInt != 1) { 129b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 13004d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 13104d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 13204d1ad83SMichael Lange } 1330877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 134*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr); 13504d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 136331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 137331830f3SMatthew G. Knepley /* Read vertices */ 138*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr); 13904d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 140331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 141*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr); 1420877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 1430877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 144854ce69bSBarry Smith ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); 145331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 146331830f3SMatthew G. Knepley int i; 147*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, &n, PETSC_ENUM);CHKERRQ(ierr); 148b9eae255SMichael Lange if (bswap) ierr = PetscByteSwap(&i, PETSC_ENUM, 1);CHKERRQ(ierr); 149*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), &m, PETSC_DOUBLE);CHKERRQ(ierr); 15004d1ad83SMichael Lange if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr); 151b9eae255SMichael Lange if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 152331830f3SMatthew G. Knepley } 153*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr);; 15404d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 155331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 156331830f3SMatthew G. Knepley /* Read cells */ 157*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr);; 15804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 159331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 160*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr);; 1610877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 1620877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 163331830f3SMatthew G. Knepley } 164db397164SMichael Lange 165abc86ac4SMichael Lange if (!rank || binary) { 166*f8e4bde8SMatthew G. Knepley PetscInt n = 1; 167a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 168a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 169a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 170a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 171a4bb7517SMichael Lange ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr); 172a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 17304d1ad83SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr); 174a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 175a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 176db397164SMichael Lange } 177*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, &n, PETSC_STRING);CHKERRQ(ierr);; 17804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 179331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 180331830f3SMatthew G. Knepley } 181abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 182abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 183a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 184331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 185331830f3SMatthew G. Knepley if (!rank) { 186a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 187a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 188a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 189a4bb7517SMichael Lange cell++; 190331830f3SMatthew G. Knepley } 191331830f3SMatthew G. Knepley } 192331830f3SMatthew G. Knepley } 193331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 194a4bb7517SMichael Lange /* Add cell-vertex connections */ 195331830f3SMatthew G. Knepley if (!rank) { 196331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 197a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 198a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 199a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 200a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1; 201331830f3SMatthew G. Knepley } 202a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 203a4bb7517SMichael Lange cell++; 204331830f3SMatthew G. Knepley } 205a4bb7517SMichael Lange } 206331830f3SMatthew G. Knepley } 2076228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 208c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 209331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 210331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 211331830f3SMatthew G. Knepley if (interpolate) { 212e56d480eSMatthew G. Knepley DM idm = NULL; 213331830f3SMatthew G. Knepley 214331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 215331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 216331830f3SMatthew G. Knepley *dm = idm; 217331830f3SMatthew G. Knepley } 21816ad7e67SMichael Lange 21916ad7e67SMichael Lange if (!rank) { 22016ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 22116ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 222d1a54cefSMatthew G. Knepley 22316ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 22416ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 225a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 22616ad7e67SMichael Lange PetscInt joinSize; 22716ad7e67SMichael Lange const PetscInt *join; 228a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 229a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1; 23016ad7e67SMichael Lange } 231a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 232a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 233a4bb7517SMichael Lange ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 234a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 23516ad7e67SMichael Lange } 236a4bb7517SMichael Lange } 23716ad7e67SMichael Lange } 23816ad7e67SMichael Lange 239331830f3SMatthew G. Knepley /* Read coordinates */ 240979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 241331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 242331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 24375b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 24475b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 245331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 246331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 247331830f3SMatthew G. Knepley } 248331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 249331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 250331830f3SMatthew G. Knepley ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 251331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 252331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 253331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 254331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 255331830f3SMatthew G. Knepley if (!rank) { 256331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 257331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 258331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 259331830f3SMatthew G. Knepley } 260331830f3SMatthew G. Knepley } 261331830f3SMatthew G. Knepley } 262331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 263331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 264331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 265331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 266a4bb7517SMichael Lange /* Clean up intermediate storage */ 267abc86ac4SMichael Lange if (!rank || binary) { 268a4bb7517SMichael Lange for (c = 0; c < numCells; ++c) { 269a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem[c].nodes); 270a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem[c].tags); 271a4bb7517SMichael Lange } 272a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem); 273a4bb7517SMichael Lange } 274331830f3SMatthew G. Knepley PetscFunctionReturn(0); 275331830f3SMatthew G. Knepley } 276db397164SMichael Lange 277db397164SMichael Lange #undef __FUNCT__ 27830412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 27904d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele) 280db397164SMichael Lange { 28104d1ad83SMichael Lange int cellType, numElem; 282*f8e4bde8SMatthew G. Knepley PetscInt n = 1; 283a4bb7517SMichael Lange PetscErrorCode ierr; 284db397164SMichael Lange 28530412aabSMichael Lange PetscFunctionBegin; 28604d1ad83SMichael Lange if (binary) { 287*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &cellType, &n, PETSC_ENUM);CHKERRQ(ierr); 288b9eae255SMichael Lange if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_ENUM, 1);CHKERRQ(ierr); 289*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &numElem, &n, PETSC_ENUM);CHKERRQ(ierr); 290b9eae255SMichael Lange if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_ENUM, 1);CHKERRQ(ierr); 291*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(ele->numTags), &n, PETSC_ENUM);CHKERRQ(ierr); 292b9eae255SMichael Lange if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_ENUM, 1);CHKERRQ(ierr); 293*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(ele->id), &n, PETSC_ENUM);CHKERRQ(ierr); 294b9eae255SMichael Lange if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_ENUM, 1);CHKERRQ(ierr); 29504d1ad83SMichael Lange } else { 296*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(ele->id), &n, PETSC_ENUM);CHKERRQ(ierr); 297*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &cellType, &n, PETSC_ENUM);CHKERRQ(ierr); 298*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(ele->numTags), &n, PETSC_ENUM);CHKERRQ(ierr); 29904d1ad83SMichael Lange } 300a4bb7517SMichael Lange ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr); 301*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, ele->tags, &ele->numTags, PETSC_ENUM);CHKERRQ(ierr); 302b9eae255SMichael Lange if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_ENUM, ele->numTags);CHKERRQ(ierr); 30330412aabSMichael Lange switch (cellType) { 30430412aabSMichael Lange case 1: /* 2-node line */ 305a4bb7517SMichael Lange ele->dim = 1; 306a4bb7517SMichael Lange ele->numNodes = 2; 30730412aabSMichael Lange break; 30830412aabSMichael Lange case 2: /* 3-node triangle */ 309a4bb7517SMichael Lange ele->dim = 2; 310a4bb7517SMichael Lange ele->numNodes = 3; 31130412aabSMichael Lange break; 31230412aabSMichael Lange case 3: /* 4-node quadrangle */ 313a4bb7517SMichael Lange ele->dim = 2; 314a4bb7517SMichael Lange ele->numNodes = 4; 31530412aabSMichael Lange break; 31630412aabSMichael Lange case 4: /* 4-node tetrahedron */ 317a4bb7517SMichael Lange ele->dim = 3; 318a4bb7517SMichael Lange ele->numNodes = 4; 31930412aabSMichael Lange break; 32030412aabSMichael Lange case 5: /* 8-node hexahedron */ 321a4bb7517SMichael Lange ele->dim = 3; 322a4bb7517SMichael Lange ele->numNodes = 8; 32330412aabSMichael Lange break; 3246b7f3382SJustin Chang case 15: /* 1-node vertex */ 325a4bb7517SMichael Lange ele->dim = 0; 326a4bb7517SMichael Lange ele->numNodes = 1; 3276b7f3382SJustin Chang break; 32830412aabSMichael Lange default: 32930412aabSMichael Lange SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 33030412aabSMichael Lange } 33104d1ad83SMichael Lange ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr); 332*f8e4bde8SMatthew G. Knepley ierr = PetscViewerRead(viewer, ele->nodes, &ele->numNodes, PETSC_ENUM);CHKERRQ(ierr); 3330877b519SMichael Lange if (byteSwap) {ierr = PetscByteSwap(ele->nodes, PETSC_ENUM, ele->numNodes);CHKERRQ(ierr);} 33430412aabSMichael Lange PetscFunctionReturn(0); 335db397164SMichael Lange } 336