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; 96*1d64f2ccSMatthew G. Knepley PetscInt dim = 0, embedDim, coordSize, c, v, d, r, cell; 97*1d64f2ccSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 98331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 99331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 100*1d64f2ccSMatthew G. Knepley PetscBool zerobase = PETSC_FALSE, isbd = PETSC_FALSE, match, binary, bswap = PETSC_FALSE; 101331830f3SMatthew G. Knepley PetscErrorCode ierr; 102331830f3SMatthew G. Knepley 103331830f3SMatthew G. Knepley PetscFunctionBegin; 104*1d64f2ccSMatthew G. Knepley ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 105*1d64f2ccSMatthew G. Knepley if (zerobase) shift = 0; 106331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 107331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 108331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 109331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1103b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 11104d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 11204d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 113abc86ac4SMichael Lange if (!rank || binary) { 114331830f3SMatthew G. Knepley PetscBool match; 115331830f3SMatthew G. Knepley int fileType, dataSize; 116f6ac7a6aSMichael Lange float version; 117331830f3SMatthew G. Knepley 118331830f3SMatthew G. Knepley /* Read header */ 119060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 12004d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 121331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 122060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 123f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 124f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 125f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 126331830f3SMatthew 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); 12704d1ad83SMichael Lange if (binary) { 128b9eae255SMichael Lange int checkInt; 129060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 13004d1ad83SMichael Lange if (checkInt != 1) { 131b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 13204d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 13304d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 13404d1ad83SMichael Lange } 1350877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 136060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 13704d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 138331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 139dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 140060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 141dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 142dc0ede3bSMatthew G. Knepley if (match) { 143dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 144dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 145dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 146dc0ede3bSMatthew G. Knepley for (r = 0; r < numRegions; ++r) { 147dc0ede3bSMatthew G. Knepley PetscInt rdim, tag; 148dc0ede3bSMatthew G. Knepley 149dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, &rdim, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 150dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, &tag, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 151dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 152dc0ede3bSMatthew G. Knepley } 153dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 154dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 155dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 156dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 157dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 158dc0ede3bSMatthew G. Knepley } 159dc0ede3bSMatthew G. Knepley /* Read vertices */ 16004d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 161331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 162060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 1630877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 1640877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 165854ce69bSBarry Smith ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); 16613aecfe5SMichael Lange if (binary) { 16713aecfe5SMichael Lange size_t doubleSize, intSize; 16813aecfe5SMichael Lange PetscInt elementSize; 16913aecfe5SMichael Lange char *buffer; 17013aecfe5SMichael Lange PetscScalar *baseptr; 17113aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); 17213aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); 17313aecfe5SMichael Lange elementSize = (intSize + 3*doubleSize); 17413aecfe5SMichael Lange ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); 17513aecfe5SMichael Lange ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); 17613aecfe5SMichael Lange if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); 177331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 17813aecfe5SMichael Lange baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize)); 17913aecfe5SMichael Lange coordsIn[v*3+0] = baseptr[0]; 18013aecfe5SMichael Lange coordsIn[v*3+1] = baseptr[1]; 18113aecfe5SMichael Lange coordsIn[v*3+2] = baseptr[2]; 18213aecfe5SMichael Lange } 18313aecfe5SMichael Lange ierr = PetscFree(buffer);CHKERRQ(ierr); 18413aecfe5SMichael Lange } else { 18513aecfe5SMichael Lange for (v = 0; v < numVertices; ++v) { 186060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 187060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 188*1d64f2ccSMatthew G. Knepley if (i != (int) v+shift) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+shift); 189331830f3SMatthew G. Knepley } 19013aecfe5SMichael Lange } 191060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 19204d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 193331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 194331830f3SMatthew G. Knepley /* Read cells */ 195060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 19604d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 197331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 198060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 1990877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 2000877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 201331830f3SMatthew G. Knepley } 202db397164SMichael Lange 203abc86ac4SMichael Lange if (!rank || binary) { 204a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 205a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 206a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 207a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 2083d51f2daSMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); 209a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 210a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 211a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 212db397164SMichael Lange } 213060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 21404d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 215331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 216331830f3SMatthew G. Knepley } 217abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 218abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 219a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 220331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 221331830f3SMatthew G. Knepley if (!rank) { 222a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 223a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 224a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 225a4bb7517SMichael Lange cell++; 226331830f3SMatthew G. Knepley } 227331830f3SMatthew G. Knepley } 228331830f3SMatthew G. Knepley } 229331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 230a4bb7517SMichael Lange /* Add cell-vertex connections */ 231331830f3SMatthew G. Knepley if (!rank) { 232331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 233a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 234a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 235a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 236*1d64f2ccSMatthew G. Knepley pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells - shift; 237331830f3SMatthew G. Knepley } 23897ed6be6Ssarens if (dim == 3) { 23997ed6be6Ssarens /* Tetrahedra are inverted */ 24097ed6be6Ssarens if (gmsh_elem[c].numNodes == 4) { 24197ed6be6Ssarens PetscInt tmp = pcone[0]; 24297ed6be6Ssarens pcone[0] = pcone[1]; 24397ed6be6Ssarens pcone[1] = tmp; 24497ed6be6Ssarens } 24597ed6be6Ssarens /* Hexahedra are inverted */ 24697ed6be6Ssarens if (gmsh_elem[c].numNodes == 8) { 24797ed6be6Ssarens PetscInt tmp = pcone[1]; 24897ed6be6Ssarens pcone[1] = pcone[3]; 24997ed6be6Ssarens pcone[3] = tmp; 25097ed6be6Ssarens } 25197ed6be6Ssarens } 252a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 253a4bb7517SMichael Lange cell++; 254331830f3SMatthew G. Knepley } 255a4bb7517SMichael Lange } 256331830f3SMatthew G. Knepley } 2576228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 258c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 259331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 260331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 261331830f3SMatthew G. Knepley if (interpolate) { 262e56d480eSMatthew G. Knepley DM idm = NULL; 263331830f3SMatthew G. Knepley 264331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 265331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 266331830f3SMatthew G. Knepley *dm = idm; 267331830f3SMatthew G. Knepley } 26816ad7e67SMichael Lange 26916ad7e67SMichael Lange if (!rank) { 27016ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 27116ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 272d1a54cefSMatthew G. Knepley 27316ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 27416ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 275a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 27616ad7e67SMichael Lange PetscInt joinSize; 27716ad7e67SMichael Lange const PetscInt *join; 278a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 279*1d64f2ccSMatthew G. Knepley pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - shift; 28016ad7e67SMichael Lange } 281a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 282a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 283c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 284a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 28516ad7e67SMichael Lange } 286a4bb7517SMichael Lange } 2876e1dd89aSLawrence Mitchell 2886e1dd89aSLawrence Mitchell /* Create cell sets */ 2896e1dd89aSLawrence Mitchell for (cell = 0, c = 0; c < numCells; ++c) { 2906e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 2916e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 2926e1dd89aSLawrence Mitchell ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 2936e1dd89aSLawrence Mitchell cell++; 2946e1dd89aSLawrence Mitchell } 2956e1dd89aSLawrence Mitchell } 2966e1dd89aSLawrence Mitchell } 29716ad7e67SMichael Lange } 29816ad7e67SMichael Lange 299331830f3SMatthew G. Knepley /* Read coordinates */ 300*1d64f2ccSMatthew G. Knepley embedDim = dim; 301*1d64f2ccSMatthew G. Knepley ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr); 302*1d64f2ccSMatthew G. Knepley if (isbd) embedDim = dim+1; 303979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 304331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 305*1d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 30675b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 30775b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 308*1d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 309*1d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 310331830f3SMatthew G. Knepley } 311331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 312331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 3138b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 314331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 315331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 316*1d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 317331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 318331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 319331830f3SMatthew G. Knepley if (!rank) { 320331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 321*1d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 322*1d64f2ccSMatthew G. Knepley coords[v*embedDim+d] = coordsIn[v*3+d]; 323331830f3SMatthew G. Knepley } 324331830f3SMatthew G. Knepley } 325331830f3SMatthew G. Knepley } 326331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 327331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 328331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 329331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 330a4bb7517SMichael Lange /* Clean up intermediate storage */ 33129403c87SMichael Lange if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 3323b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 333331830f3SMatthew G. Knepley PetscFunctionReturn(0); 334331830f3SMatthew G. Knepley } 335db397164SMichael Lange 336db397164SMichael Lange #undef __FUNCT__ 33730412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 3383d51f2daSMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems) 339db397164SMichael Lange { 3405257f852SMichael Lange PetscInt c, p; 3413d51f2daSMichael Lange GmshElement *elements; 342b6a3fe3cSMichael Lange int i, cellType, dim, numNodes, numElem, numTags; 3435257f852SMichael Lange int ibuf[16]; 344a4bb7517SMichael Lange PetscErrorCode ierr; 345db397164SMichael Lange 34630412aabSMichael Lange PetscFunctionBegin; 3473d51f2daSMichael Lange ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 3483d51f2daSMichael Lange for (c = 0; c < numCells;) { 349b6a3fe3cSMichael Lange ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 350b6a3fe3cSMichael Lange if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr); 35104d1ad83SMichael Lange if (binary) { 352b6a3fe3cSMichael Lange cellType = ibuf[0]; 353b6a3fe3cSMichael Lange numElem = ibuf[1]; 354b6a3fe3cSMichael Lange numTags = ibuf[2]; 35504d1ad83SMichael Lange } else { 356b6a3fe3cSMichael Lange elements[c].id = ibuf[0]; 357b6a3fe3cSMichael Lange cellType = ibuf[1]; 358b6a3fe3cSMichael Lange numTags = ibuf[2]; 3593d51f2daSMichael Lange numElem = 1; 36004d1ad83SMichael Lange } 361b6a3fe3cSMichael Lange switch (cellType) { 362b6a3fe3cSMichael Lange case 1: /* 2-node line */ 363b6a3fe3cSMichael Lange dim = 1; 364b6a3fe3cSMichael Lange numNodes = 2; 365b6a3fe3cSMichael Lange break; 366b6a3fe3cSMichael Lange case 2: /* 3-node triangle */ 367b6a3fe3cSMichael Lange dim = 2; 368b6a3fe3cSMichael Lange numNodes = 3; 369b6a3fe3cSMichael Lange break; 370b6a3fe3cSMichael Lange case 3: /* 4-node quadrangle */ 371b6a3fe3cSMichael Lange dim = 2; 372b6a3fe3cSMichael Lange numNodes = 4; 373b6a3fe3cSMichael Lange break; 374b6a3fe3cSMichael Lange case 4: /* 4-node tetrahedron */ 375b6a3fe3cSMichael Lange dim = 3; 376b6a3fe3cSMichael Lange numNodes = 4; 377b6a3fe3cSMichael Lange break; 378b6a3fe3cSMichael Lange case 5: /* 8-node hexahedron */ 379b6a3fe3cSMichael Lange dim = 3; 380b6a3fe3cSMichael Lange numNodes = 8; 381b6a3fe3cSMichael Lange break; 382b6a3fe3cSMichael Lange case 15: /* 1-node vertex */ 383b6a3fe3cSMichael Lange dim = 0; 384b6a3fe3cSMichael Lange numNodes = 1; 385b6a3fe3cSMichael Lange break; 386b6a3fe3cSMichael Lange default: 387b6a3fe3cSMichael Lange SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 388b6a3fe3cSMichael Lange } 3895257f852SMichael Lange if (binary) { 3905257f852SMichael Lange const PetscInt nint = numNodes + numTags + 1; 3913d51f2daSMichael Lange for (i = 0; i < numElem; ++i, ++c) { 3923d51f2daSMichael Lange /* Loop over inner binary element block */ 393b6a3fe3cSMichael Lange elements[c].dim = dim; 394b6a3fe3cSMichael Lange elements[c].numNodes = numNodes; 395b6a3fe3cSMichael Lange elements[c].numTags = numTags; 396b6a3fe3cSMichael Lange 3975257f852SMichael Lange ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 3985257f852SMichael Lange if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr); 3995257f852SMichael Lange elements[c].id = ibuf[0]; 4005257f852SMichael Lange for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 4015257f852SMichael Lange for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 4025257f852SMichael Lange } 4035257f852SMichael Lange } else { 4045257f852SMichael Lange elements[c].dim = dim; 4055257f852SMichael Lange elements[c].numNodes = numNodes; 4065257f852SMichael Lange elements[c].numTags = numTags; 4073d51f2daSMichael Lange ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 4083d51f2daSMichael Lange ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 4095257f852SMichael Lange c++; 4103d51f2daSMichael Lange } 4113d51f2daSMichael Lange } 4123d51f2daSMichael Lange *gmsh_elems = elements; 41330412aabSMichael Lange PetscFunctionReturn(0); 414db397164SMichael Lange } 415