1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #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) { 417d282ae0SMichael Lange /* Read only the first two lines of the Gmsh file */ 42060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 437d282ae0SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 447d282ae0SMichael Lange if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 45060da220SMatthew G. Knepley ierr = PetscViewerRead(vheader, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 46c1c22fd2SMatthew G. Knepley snum = sscanf(line, "%f %d", &version, &fT); 47c1c22fd2SMatthew G. Knepley fileType = (PetscInt) fT; 48f6ac7a6aSMichael Lange if (snum != 2) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 49f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 50abc86ac4SMichael Lange } 51abc86ac4SMichael Lange ierr = MPI_Bcast(&fileType, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 527d282ae0SMichael Lange /* Create appropriate viewer and build plex */ 537d282ae0SMichael Lange if (fileType == 0) vtype = PETSCVIEWERASCII; 547d282ae0SMichael Lange else vtype = PETSCVIEWERBINARY; 557d282ae0SMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 567d282ae0SMichael Lange ierr = PetscViewerSetType(viewer, vtype);CHKERRQ(ierr); 577d282ae0SMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 587d282ae0SMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 597d282ae0SMichael Lange ierr = DMPlexCreateGmsh(comm, viewer, interpolate, dm);CHKERRQ(ierr); 607d282ae0SMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 617d282ae0SMichael Lange ierr = PetscViewerDestroy(&vheader);CHKERRQ(ierr); 627d282ae0SMichael Lange PetscFunctionReturn(0); 637d282ae0SMichael Lange } 647d282ae0SMichael Lange 657d282ae0SMichael Lange 667d282ae0SMichael Lange #undef __FUNCT__ 67331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh" 68331830f3SMatthew G. Knepley /*@ 697d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 70331830f3SMatthew G. Knepley 71331830f3SMatthew G. Knepley Collective on comm 72331830f3SMatthew G. Knepley 73331830f3SMatthew G. Knepley Input Parameters: 74331830f3SMatthew G. Knepley + comm - The MPI communicator 75331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 76331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 77331830f3SMatthew G. Knepley 78331830f3SMatthew G. Knepley Output Parameter: 79331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 80331830f3SMatthew G. Knepley 81331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 823d51f2daSMichael Lange and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-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; 9713aecfe5SMichael Lange int i, 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); 1083b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 10904d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 11004d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 111abc86ac4SMichael Lange if (!rank || binary) { 112331830f3SMatthew G. Knepley PetscBool match; 113331830f3SMatthew G. Knepley int fileType, dataSize; 114f6ac7a6aSMichael Lange float version; 115331830f3SMatthew G. Knepley 116331830f3SMatthew G. Knepley /* Read header */ 117060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, 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"); 120060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, 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; 127060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, 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); 134060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, 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 */ 138060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, 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"); 141060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, 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); 14513aecfe5SMichael Lange if (binary) { 14613aecfe5SMichael Lange size_t doubleSize, intSize; 14713aecfe5SMichael Lange PetscInt elementSize; 14813aecfe5SMichael Lange char *buffer; 14913aecfe5SMichael Lange PetscScalar *baseptr; 15013aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); 15113aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); 15213aecfe5SMichael Lange elementSize = (intSize + 3*doubleSize); 15313aecfe5SMichael Lange ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); 15413aecfe5SMichael Lange ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); 15513aecfe5SMichael Lange if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); 156331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 15713aecfe5SMichael Lange baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize)); 15813aecfe5SMichael Lange coordsIn[v*3+0] = baseptr[0]; 15913aecfe5SMichael Lange coordsIn[v*3+1] = baseptr[1]; 16013aecfe5SMichael Lange coordsIn[v*3+2] = baseptr[2]; 16113aecfe5SMichael Lange } 16213aecfe5SMichael Lange ierr = PetscFree(buffer);CHKERRQ(ierr); 16313aecfe5SMichael Lange } else { 16413aecfe5SMichael Lange for (v = 0; v < numVertices; ++v) { 165060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 166060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 167b9eae255SMichael Lange if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 168331830f3SMatthew G. Knepley } 16913aecfe5SMichael Lange } 170060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 17104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 172331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 173331830f3SMatthew G. Knepley /* Read cells */ 174060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 17504d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 176331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 177060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 1780877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 1790877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 180331830f3SMatthew G. Knepley } 181db397164SMichael Lange 182abc86ac4SMichael Lange if (!rank || binary) { 183a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 184a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 185a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 186a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 1873d51f2daSMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); 188a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 189a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 190a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 191db397164SMichael Lange } 192060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 19304d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 194331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 195331830f3SMatthew G. Knepley } 196abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 197abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 198a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 199331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 200331830f3SMatthew G. Knepley if (!rank) { 201a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 202a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 203a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 204a4bb7517SMichael Lange cell++; 205331830f3SMatthew G. Knepley } 206331830f3SMatthew G. Knepley } 207331830f3SMatthew G. Knepley } 208331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 209a4bb7517SMichael Lange /* Add cell-vertex connections */ 210331830f3SMatthew G. Knepley if (!rank) { 211331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 212a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 213a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 214a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 215a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1; 216331830f3SMatthew G. Knepley } 217a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 218a4bb7517SMichael Lange cell++; 219331830f3SMatthew G. Knepley } 220a4bb7517SMichael Lange } 221331830f3SMatthew G. Knepley } 2226228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 223c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 224331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 225331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 226331830f3SMatthew G. Knepley if (interpolate) { 227e56d480eSMatthew G. Knepley DM idm = NULL; 228331830f3SMatthew G. Knepley 229331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 230331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 231331830f3SMatthew G. Knepley *dm = idm; 232331830f3SMatthew G. Knepley } 23316ad7e67SMichael Lange 23416ad7e67SMichael Lange if (!rank) { 23516ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 23616ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 237d1a54cefSMatthew G. Knepley 23816ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 23916ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 240a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 24116ad7e67SMichael Lange PetscInt joinSize; 24216ad7e67SMichael Lange const PetscInt *join; 243a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 244a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1; 24516ad7e67SMichael Lange } 246a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 247a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 248c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 249a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 25016ad7e67SMichael Lange } 251a4bb7517SMichael Lange } 25216ad7e67SMichael Lange } 25316ad7e67SMichael Lange 254331830f3SMatthew G. Knepley /* Read coordinates */ 255979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 256331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 257331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 25875b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 25975b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 260331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 261331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 262331830f3SMatthew G. Knepley } 263331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 264331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 265*8b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 266331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 267331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 268331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 269331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 270331830f3SMatthew G. Knepley if (!rank) { 271331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 272331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 273331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 274331830f3SMatthew G. Knepley } 275331830f3SMatthew G. Knepley } 276331830f3SMatthew G. Knepley } 277331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 278331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 279331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 280331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 281a4bb7517SMichael Lange /* Clean up intermediate storage */ 28229403c87SMichael Lange if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 2833b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 284331830f3SMatthew G. Knepley PetscFunctionReturn(0); 285331830f3SMatthew G. Knepley } 286db397164SMichael Lange 287db397164SMichael Lange #undef __FUNCT__ 28830412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 2893d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems) 290db397164SMichael Lange { 2915257f852SMichael Lange PetscInt c, p; 2923d51f2daSMichael Lange GmshElement *elements; 293b6a3fe3cSMichael Lange int i, cellType, dim, numNodes, numElem, numTags; 2945257f852SMichael Lange int ibuf[16]; 295a4bb7517SMichael Lange PetscErrorCode ierr; 296db397164SMichael Lange 29730412aabSMichael Lange PetscFunctionBegin; 2983d51f2daSMichael Lange ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 2993d51f2daSMichael Lange for (c = 0; c < numCells;) { 300b6a3fe3cSMichael Lange ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 301b6a3fe3cSMichael Lange if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr); 30204d1ad83SMichael Lange if (binary) { 303b6a3fe3cSMichael Lange cellType = ibuf[0]; 304b6a3fe3cSMichael Lange numElem = ibuf[1]; 305b6a3fe3cSMichael Lange numTags = ibuf[2]; 30604d1ad83SMichael Lange } else { 307b6a3fe3cSMichael Lange elements[c].id = ibuf[0]; 308b6a3fe3cSMichael Lange cellType = ibuf[1]; 309b6a3fe3cSMichael Lange numTags = ibuf[2]; 3103d51f2daSMichael Lange numElem = 1; 31104d1ad83SMichael Lange } 312b6a3fe3cSMichael Lange switch (cellType) { 313b6a3fe3cSMichael Lange case 1: /* 2-node line */ 314b6a3fe3cSMichael Lange dim = 1; 315b6a3fe3cSMichael Lange numNodes = 2; 316b6a3fe3cSMichael Lange break; 317b6a3fe3cSMichael Lange case 2: /* 3-node triangle */ 318b6a3fe3cSMichael Lange dim = 2; 319b6a3fe3cSMichael Lange numNodes = 3; 320b6a3fe3cSMichael Lange break; 321b6a3fe3cSMichael Lange case 3: /* 4-node quadrangle */ 322b6a3fe3cSMichael Lange dim = 2; 323b6a3fe3cSMichael Lange numNodes = 4; 324b6a3fe3cSMichael Lange break; 325b6a3fe3cSMichael Lange case 4: /* 4-node tetrahedron */ 326b6a3fe3cSMichael Lange dim = 3; 327b6a3fe3cSMichael Lange numNodes = 4; 328b6a3fe3cSMichael Lange break; 329b6a3fe3cSMichael Lange case 5: /* 8-node hexahedron */ 330b6a3fe3cSMichael Lange dim = 3; 331b6a3fe3cSMichael Lange numNodes = 8; 332b6a3fe3cSMichael Lange break; 333b6a3fe3cSMichael Lange case 15: /* 1-node vertex */ 334b6a3fe3cSMichael Lange dim = 0; 335b6a3fe3cSMichael Lange numNodes = 1; 336b6a3fe3cSMichael Lange break; 337b6a3fe3cSMichael Lange default: 338b6a3fe3cSMichael Lange SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 339b6a3fe3cSMichael Lange } 3405257f852SMichael Lange if (binary) { 3415257f852SMichael Lange const PetscInt nint = numNodes + numTags + 1; 3423d51f2daSMichael Lange for (i = 0; i < numElem; ++i, ++c) { 3433d51f2daSMichael Lange /* Loop over inner binary element block */ 344b6a3fe3cSMichael Lange elements[c].dim = dim; 345b6a3fe3cSMichael Lange elements[c].numNodes = numNodes; 346b6a3fe3cSMichael Lange elements[c].numTags = numTags; 347b6a3fe3cSMichael Lange 3485257f852SMichael Lange ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 3495257f852SMichael Lange if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr); 3505257f852SMichael Lange elements[c].id = ibuf[0]; 3515257f852SMichael Lange for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 3525257f852SMichael Lange for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 3535257f852SMichael Lange } 3545257f852SMichael Lange } else { 3555257f852SMichael Lange elements[c].dim = dim; 3565257f852SMichael Lange elements[c].numNodes = numNodes; 3575257f852SMichael Lange elements[c].numTags = numTags; 3583d51f2daSMichael Lange ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 3593d51f2daSMichael Lange ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 3605257f852SMichael Lange c++; 3613d51f2daSMichael Lange } 3623d51f2daSMichael Lange } 3633d51f2daSMichael Lange *gmsh_elems = elements; 36430412aabSMichael Lange PetscFunctionReturn(0); 365db397164SMichael Lange } 366