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 63df032b82SMatthew G. Knepley static PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscInt numCells, PetscBool binary, PetscBool byteSwap, GmshElement **gmsh_elems) 64df032b82SMatthew G. Knepley { 65df032b82SMatthew G. Knepley PetscInt c, p; 66df032b82SMatthew G. Knepley GmshElement *elements; 67bf6ba3a3SMatthew G. Knepley int i, cellType, dim, numNodes, numNodesIgnore, numElem, numTags; 68bf6ba3a3SMatthew G. Knepley PetscInt pibuf[64]; 69df032b82SMatthew G. Knepley int ibuf[16]; 70df032b82SMatthew G. Knepley PetscErrorCode ierr; 71df032b82SMatthew G. Knepley 72df032b82SMatthew G. Knepley PetscFunctionBegin; 73df032b82SMatthew G. Knepley ierr = PetscMalloc1(numCells, &elements);CHKERRQ(ierr); 74df032b82SMatthew G. Knepley for (c = 0; c < numCells;) { 75df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, &ibuf, 3, NULL, PETSC_ENUM);CHKERRQ(ierr); 76df032b82SMatthew G. Knepley if (byteSwap) ierr = PetscByteSwap(&ibuf, PETSC_ENUM, 3);CHKERRQ(ierr); 77df032b82SMatthew G. Knepley if (binary) { 78df032b82SMatthew G. Knepley cellType = ibuf[0]; 79df032b82SMatthew G. Knepley numElem = ibuf[1]; 80df032b82SMatthew G. Knepley numTags = ibuf[2]; 81df032b82SMatthew G. Knepley } else { 82df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 83df032b82SMatthew G. Knepley cellType = ibuf[1]; 84df032b82SMatthew G. Knepley numTags = ibuf[2]; 85df032b82SMatthew G. Knepley numElem = 1; 86df032b82SMatthew G. Knepley } 87bf6ba3a3SMatthew G. Knepley /* http://gmsh.info/doc/texinfo/gmsh.html#MSH-ASCII-file-format */ 88bf6ba3a3SMatthew G. Knepley numNodesIgnore = 0; 89df032b82SMatthew G. Knepley switch (cellType) { 90df032b82SMatthew G. Knepley case 1: /* 2-node line */ 91df032b82SMatthew G. Knepley dim = 1; 92df032b82SMatthew G. Knepley numNodes = 2; 93df032b82SMatthew G. Knepley break; 94df032b82SMatthew G. Knepley case 2: /* 3-node triangle */ 95df032b82SMatthew G. Knepley dim = 2; 96df032b82SMatthew G. Knepley numNodes = 3; 97df032b82SMatthew G. Knepley break; 98df032b82SMatthew G. Knepley case 3: /* 4-node quadrangle */ 99df032b82SMatthew G. Knepley dim = 2; 100df032b82SMatthew G. Knepley numNodes = 4; 101df032b82SMatthew G. Knepley break; 102df032b82SMatthew G. Knepley case 4: /* 4-node tetrahedron */ 103df032b82SMatthew G. Knepley dim = 3; 104df032b82SMatthew G. Knepley numNodes = 4; 105df032b82SMatthew G. Knepley break; 106df032b82SMatthew G. Knepley case 5: /* 8-node hexahedron */ 107df032b82SMatthew G. Knepley dim = 3; 108df032b82SMatthew G. Knepley numNodes = 8; 109df032b82SMatthew G. Knepley break; 110bf6ba3a3SMatthew G. Knepley case 8: /* 3-node 2nd order line */ 111bf6ba3a3SMatthew G. Knepley dim = 1; 112bf6ba3a3SMatthew G. Knepley numNodes = 2; 113bf6ba3a3SMatthew G. Knepley numNodesIgnore = 1; 114bf6ba3a3SMatthew G. Knepley break; 115bf6ba3a3SMatthew G. Knepley case 9: /* 6-node 2nd order triangle */ 116bf6ba3a3SMatthew G. Knepley dim = 2; 117bf6ba3a3SMatthew G. Knepley numNodes = 3; 118bf6ba3a3SMatthew G. Knepley numNodesIgnore = 3; 119bf6ba3a3SMatthew G. Knepley break; 120df032b82SMatthew G. Knepley case 15: /* 1-node vertex */ 121df032b82SMatthew G. Knepley dim = 0; 122df032b82SMatthew G. Knepley numNodes = 1; 123df032b82SMatthew G. Knepley break; 124bf6ba3a3SMatthew G. Knepley case 6: /* 6-node prism */ 125bf6ba3a3SMatthew G. Knepley case 7: /* 5-node pyramid */ 126bf6ba3a3SMatthew G. Knepley case 10: /* 9-node 2nd order quadrangle */ 127bf6ba3a3SMatthew G. Knepley case 11: /* 10-node 2nd order tetrahedron */ 128bf6ba3a3SMatthew G. Knepley case 12: /* 27-node 2nd order hexhedron */ 129bf6ba3a3SMatthew G. Knepley case 13: /* 19-node 2nd order prism */ 130bf6ba3a3SMatthew G. Knepley case 14: /* 14-node 2nd order pyramid */ 131df032b82SMatthew G. Knepley default: 132df032b82SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 133df032b82SMatthew G. Knepley } 134df032b82SMatthew G. Knepley if (binary) { 135bf6ba3a3SMatthew G. Knepley const PetscInt nint = numNodes + numTags + 1 + numNodesIgnore; 136df032b82SMatthew G. Knepley for (i = 0; i < numElem; ++i, ++c) { 137df032b82SMatthew G. Knepley /* Loop over inner binary element block */ 138df032b82SMatthew G. Knepley elements[c].dim = dim; 139df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 140df032b82SMatthew G. Knepley elements[c].numTags = numTags; 141df032b82SMatthew G. Knepley 142df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, &ibuf, nint, NULL, PETSC_ENUM);CHKERRQ(ierr); 143df032b82SMatthew G. Knepley if (byteSwap) ierr = PetscByteSwap( &ibuf, PETSC_ENUM, nint);CHKERRQ(ierr); 144df032b82SMatthew G. Knepley elements[c].id = ibuf[0]; 145df032b82SMatthew G. Knepley for (p = 0; p < numTags; p++) elements[c].tags[p] = ibuf[1 + p]; 146df032b82SMatthew G. Knepley for (p = 0; p < numNodes; p++) elements[c].nodes[p] = ibuf[1 + numTags + p]; 147df032b82SMatthew G. Knepley } 148df032b82SMatthew G. Knepley } else { 149df032b82SMatthew G. Knepley elements[c].dim = dim; 150df032b82SMatthew G. Knepley elements[c].numNodes = numNodes; 151df032b82SMatthew G. Knepley elements[c].numTags = numTags; 152df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].tags, elements[c].numTags, NULL, PETSC_ENUM);CHKERRQ(ierr); 153df032b82SMatthew G. Knepley ierr = PetscViewerRead(viewer, elements[c].nodes, elements[c].numNodes, NULL, PETSC_ENUM);CHKERRQ(ierr); 154bf6ba3a3SMatthew G. Knepley ierr = PetscViewerRead(viewer, pibuf, numNodesIgnore, NULL, PETSC_ENUM);CHKERRQ(ierr); 155df032b82SMatthew G. Knepley c++; 156df032b82SMatthew G. Knepley } 157df032b82SMatthew G. Knepley } 158df032b82SMatthew G. Knepley *gmsh_elems = elements; 159df032b82SMatthew G. Knepley PetscFunctionReturn(0); 160df032b82SMatthew G. Knepley } 1617d282ae0SMichael Lange 162331830f3SMatthew G. Knepley /*@ 1637d282ae0SMichael Lange DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file viewer 164331830f3SMatthew G. Knepley 165331830f3SMatthew G. Knepley Collective on comm 166331830f3SMatthew G. Knepley 167331830f3SMatthew G. Knepley Input Parameters: 168331830f3SMatthew G. Knepley + comm - The MPI communicator 169331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 170331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 171331830f3SMatthew G. Knepley 172331830f3SMatthew G. Knepley Output Parameter: 173331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 174331830f3SMatthew G. Knepley 175331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 1763d51f2daSMichael Lange and http://www.geuz.org/gmsh/doc/texinfo/#MSH-binary-file-format 177331830f3SMatthew G. Knepley 178331830f3SMatthew G. Knepley Level: beginner 179331830f3SMatthew G. Knepley 180331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 181331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 182331830f3SMatthew G. Knepley @*/ 183331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 184331830f3SMatthew G. Knepley { 18504d1ad83SMichael Lange PetscViewerType vtype; 1863c67d5adSSatish Balay GmshElement *gmsh_elem = NULL; 187331830f3SMatthew G. Knepley PetscSection coordSection; 188331830f3SMatthew G. Knepley Vec coordinates; 1896fbe17bfSStefano Zampini PetscBT periodicV = NULL, periodicC = NULL; 19084572febSToby Isaac PetscScalar *coords; 19184572febSToby Isaac PetscReal *coordsIn = NULL; 192*dd169d64SMatthew G. Knepley PetscInt dim = 0, embedDim, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL; 1931d64f2ccSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum, shift = 1; 194331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 195331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 196*dd169d64SMatthew G. Knepley PetscBool zerobase = PETSC_FALSE, isbd = PETSC_FALSE, match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 197331830f3SMatthew G. Knepley PetscErrorCode ierr; 198331830f3SMatthew G. Knepley 199331830f3SMatthew G. Knepley PetscFunctionBegin; 2001d64f2ccSMatthew G. Knepley ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_zero_base", &zerobase, NULL);CHKERRQ(ierr); 2011d64f2ccSMatthew G. Knepley if (zerobase) shift = 0; 202331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 203331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 204331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 205331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2063b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 20704d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 20804d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 209d3f73514SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_periodic",&periodic,NULL);CHKERRQ(ierr); 210d3f73514SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_usemarker",&usemarker,NULL);CHKERRQ(ierr); 211abc86ac4SMichael Lange if (!rank || binary) { 212331830f3SMatthew G. Knepley PetscBool match; 213331830f3SMatthew G. Knepley int fileType, dataSize; 214f6ac7a6aSMichael Lange float version; 215331830f3SMatthew G. Knepley 216331830f3SMatthew G. Knepley /* Read header */ 217060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 21804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 219331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 220060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 221f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 222f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 223f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 224331830f3SMatthew 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); 22504d1ad83SMichael Lange if (binary) { 226b9eae255SMichael Lange int checkInt; 227060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 22804d1ad83SMichael Lange if (checkInt != 1) { 229b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 23004d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 23104d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 23204d1ad83SMichael Lange } 2330877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 234060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 23504d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 236331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 237dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 238060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 239dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 240dc0ede3bSMatthew G. Knepley if (match) { 241dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 242dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 243dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 244a8667449SSander Arens for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);} 245dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 246dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 247dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 248dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 249dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 250dc0ede3bSMatthew G. Knepley } 251dc0ede3bSMatthew G. Knepley /* Read vertices */ 25204d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 253331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 254060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 2550877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 2560877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 257854ce69bSBarry Smith ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); 25813aecfe5SMichael Lange if (binary) { 25913aecfe5SMichael Lange size_t doubleSize, intSize; 26013aecfe5SMichael Lange PetscInt elementSize; 26113aecfe5SMichael Lange char *buffer; 262e15165f2SToby Isaac double *baseptr; 26313aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); 26413aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); 26513aecfe5SMichael Lange elementSize = (intSize + 3*doubleSize); 26613aecfe5SMichael Lange ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); 26713aecfe5SMichael Lange ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); 26813aecfe5SMichael Lange if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); 269331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 270e15165f2SToby Isaac baseptr = ((double*)(buffer+v*elementSize+intSize)); 271e15165f2SToby Isaac coordsIn[v*3+0] = (PetscReal) baseptr[0]; 272e15165f2SToby Isaac coordsIn[v*3+1] = (PetscReal) baseptr[1]; 273e15165f2SToby Isaac coordsIn[v*3+2] = (PetscReal) baseptr[2]; 27413aecfe5SMichael Lange } 27513aecfe5SMichael Lange ierr = PetscFree(buffer);CHKERRQ(ierr); 27613aecfe5SMichael Lange } else { 27713aecfe5SMichael Lange for (v = 0; v < numVertices; ++v) { 278060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 279060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 2801d64f2ccSMatthew 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); 281331830f3SMatthew G. Knepley } 28213aecfe5SMichael Lange } 283060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 28404d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 285331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 286331830f3SMatthew G. Knepley /* Read cells */ 287060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 28804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 289331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 290060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 2910877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 2920877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 293331830f3SMatthew G. Knepley } 294db397164SMichael Lange 295abc86ac4SMichael Lange if (!rank || binary) { 296a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 297a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 298a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 299a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 3003d51f2daSMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); 301a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 302a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 303a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 304db397164SMichael Lange } 305d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 3061b82c3ebSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 307331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 308d08df55aSStefano Zampini if (periodic) { 309f45c9589SStefano Zampini PetscInt pVert, *periodicMapT, *aux; 310d08df55aSStefano Zampini int numPeriodic; 311d08df55aSStefano Zampini 312d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 313d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 314d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 315f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 3166fbe17bfSStefano Zampini ierr = PetscBTCreate(numVertices, &periodicV);CHKERRQ(ierr); 317f45c9589SStefano Zampini for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 318d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 319d08df55aSStefano Zampini snum = sscanf(line, "%d", &numPeriodic); 320d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 321d08df55aSStefano Zampini for (i = 0; i < numPeriodic; i++) { 322d08df55aSStefano Zampini int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes; 323d08df55aSStefano Zampini 324d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 325d08df55aSStefano Zampini snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */ 326d08df55aSStefano Zampini if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 32772ffbcc9SStefano Zampini /* Type of tranformation, discarded */ 32872ffbcc9SStefano Zampini ierr = PetscViewerRead(viewer, line, -PETSC_MAX_PATH_LEN, NULL, PETSC_STRING);CHKERRQ(ierr); 329d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 330d08df55aSStefano Zampini snum = sscanf(line, "%d", &nNodes); 331d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 332d08df55aSStefano Zampini for (j = 0; j < nNodes; j++) { 333d08df55aSStefano Zampini PetscInt slaveNode, masterNode; 334d08df55aSStefano Zampini 335d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 336a942bdb0SSatish Balay snum = sscanf(line, "%" PetscInt_FMT " %" PetscInt_FMT, &slaveNode, &masterNode); 337d08df55aSStefano Zampini if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 338f45c9589SStefano Zampini periodicMapT[slaveNode - 1] = masterNode - 1; 339437bdd5bSStefano Zampini ierr = PetscBTSet(periodicV, slaveNode - 1);CHKERRQ(ierr); 340437bdd5bSStefano Zampini ierr = PetscBTSet(periodicV, masterNode - 1);CHKERRQ(ierr); 341d08df55aSStefano Zampini } 342d08df55aSStefano Zampini } 343d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 344d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 345d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 346d08df55aSStefano Zampini /* we may have slaves of slaves */ 347d08df55aSStefano Zampini for (i = 0; i < numVertices; i++) { 348f45c9589SStefano Zampini while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 349f45c9589SStefano Zampini periodicMapT[i] = periodicMapT[periodicMapT[i]]; 350d08df55aSStefano Zampini } 351d08df55aSStefano Zampini } 352f45c9589SStefano Zampini /* periodicMap : from old to new numbering (periodic vertices excluded) 353f45c9589SStefano Zampini periodicMapI: from new to old numbering */ 354f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 355f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 356f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 357f45c9589SStefano Zampini for (i = 0, pVert = 0; i < numVertices; i++) { 358f45c9589SStefano Zampini if (periodicMapT[i] != i) { 359f45c9589SStefano Zampini pVert++; 360f45c9589SStefano Zampini } else { 361f45c9589SStefano Zampini aux[i] = i - pVert; 362f45c9589SStefano Zampini periodicMapI[i - pVert] = i; 363f45c9589SStefano Zampini } 364f45c9589SStefano Zampini } 365f45c9589SStefano Zampini for (i = 0 ; i < numVertices; i++) { 366f45c9589SStefano Zampini periodicMap[i] = aux[periodicMapT[i]]; 367f45c9589SStefano Zampini } 368f45c9589SStefano Zampini ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 369f45c9589SStefano Zampini ierr = PetscFree(aux);CHKERRQ(ierr); 370f45c9589SStefano Zampini /* remove periodic vertices */ 371f45c9589SStefano Zampini numVertices = numVertices - pVert; 372d08df55aSStefano Zampini } 373331830f3SMatthew G. Knepley } 374abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 375abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 376a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 377331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 378331830f3SMatthew G. Knepley if (!rank) { 379a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 380a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 381a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 382a4bb7517SMichael Lange cell++; 383331830f3SMatthew G. Knepley } 384331830f3SMatthew G. Knepley } 385331830f3SMatthew G. Knepley } 386331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 387a4bb7517SMichael Lange /* Add cell-vertex connections */ 388331830f3SMatthew G. Knepley if (!rank) { 389331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 390437bdd5bSStefano Zampini 391a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 392a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 393a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 394*dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 3956fbe17bfSStefano Zampini pcone[corner] = periodicMap ? periodicMap[cc] + trueNumCells 3966fbe17bfSStefano Zampini : cc + trueNumCells; 397331830f3SMatthew G. Knepley } 39897ed6be6Ssarens if (dim == 3) { 39997ed6be6Ssarens /* Tetrahedra are inverted */ 40097ed6be6Ssarens if (gmsh_elem[c].numNodes == 4) { 40197ed6be6Ssarens PetscInt tmp = pcone[0]; 40297ed6be6Ssarens pcone[0] = pcone[1]; 40397ed6be6Ssarens pcone[1] = tmp; 40497ed6be6Ssarens } 40597ed6be6Ssarens /* Hexahedra are inverted */ 40697ed6be6Ssarens if (gmsh_elem[c].numNodes == 8) { 40797ed6be6Ssarens PetscInt tmp = pcone[1]; 40897ed6be6Ssarens pcone[1] = pcone[3]; 40997ed6be6Ssarens pcone[3] = tmp; 41097ed6be6Ssarens } 41197ed6be6Ssarens } 412a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 413a4bb7517SMichael Lange cell++; 414331830f3SMatthew G. Knepley } 415a4bb7517SMichael Lange } 416331830f3SMatthew G. Knepley } 4176228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 418c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 419331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 420331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 421331830f3SMatthew G. Knepley if (interpolate) { 422e56d480eSMatthew G. Knepley DM idm = NULL; 423331830f3SMatthew G. Knepley 424331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 425331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 426331830f3SMatthew G. Knepley *dm = idm; 427331830f3SMatthew G. Knepley } 428d3f73514SStefano Zampini if (usemarker) { 429d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 430d3f73514SStefano Zampini 431d3f73514SStefano Zampini ierr = DMCreateLabel(*dm,"marker");CHKERRQ(ierr); 432d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 433d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 434d3f73514SStefano Zampini PetscInt suppSize; 435d3f73514SStefano Zampini 436d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 437d3f73514SStefano Zampini if (suppSize == 1) { 438d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 439d3f73514SStefano Zampini 440d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 441d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 442d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 443d3f73514SStefano Zampini } 444d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 445d3f73514SStefano Zampini } 446d3f73514SStefano Zampini } 447d3f73514SStefano Zampini } 44816ad7e67SMichael Lange 44916ad7e67SMichael Lange if (!rank) { 45016ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 45116ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 452d1a54cefSMatthew G. Knepley 45316ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 45416ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 455a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 45616ad7e67SMichael Lange const PetscInt *join; 457437bdd5bSStefano Zampini PetscInt joinSize; 458437bdd5bSStefano Zampini 459a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 460*dd169d64SMatthew G. Knepley const PetscInt cc = gmsh_elem[c].nodes[corner] - shift; 461437bdd5bSStefano Zampini pcone[corner] = periodicMap ? periodicMap[cc] + vStart 462437bdd5bSStefano Zampini : cc + vStart; 46316ad7e67SMichael Lange } 464a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 465a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 466c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 467a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 46816ad7e67SMichael Lange } 469a4bb7517SMichael Lange } 4706e1dd89aSLawrence Mitchell 4716e1dd89aSLawrence Mitchell /* Create cell sets */ 4726e1dd89aSLawrence Mitchell for (cell = 0, c = 0; c < numCells; ++c) { 4736e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 4746e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 4756e1dd89aSLawrence Mitchell ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 4766e1dd89aSLawrence Mitchell cell++; 4776e1dd89aSLawrence Mitchell } 4786e1dd89aSLawrence Mitchell } 4796e1dd89aSLawrence Mitchell } 4800c070f12SSander Arens 4810c070f12SSander Arens /* Create vertex sets */ 4820c070f12SSander Arens for (c = 0; c < numCells; ++c) { 4830c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 4840c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 485d08df55aSStefano Zampini PetscInt vid = periodicMap ? periodicMap[gmsh_elem[c].nodes[0] - 1] + vStart 486d08df55aSStefano Zampini : gmsh_elem[c].nodes[0] - 1 + vStart; 487d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 4880c070f12SSander Arens } 4890c070f12SSander Arens } 4900c070f12SSander Arens } 49116ad7e67SMichael Lange } 49216ad7e67SMichael Lange 493331830f3SMatthew G. Knepley /* Read coordinates */ 4941d64f2ccSMatthew G. Knepley embedDim = dim; 4951d64f2ccSMatthew G. Knepley ierr = PetscOptionsGetBool(NULL, NULL, "-gmsh_bd", &isbd, NULL);CHKERRQ(ierr); 4961d64f2ccSMatthew G. Knepley if (isbd) embedDim = dim+1; 497979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 498331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 4991d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, embedDim);CHKERRQ(ierr); 500f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 501f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 502f45c9589SStefano Zampini } else { 50375b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 504f45c9589SStefano Zampini } 50575b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 5061d64f2ccSMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, embedDim);CHKERRQ(ierr); 5071d64f2ccSMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, embedDim);CHKERRQ(ierr); 508331830f3SMatthew G. Knepley } 509f45c9589SStefano Zampini if (periodicMap) { 510437bdd5bSStefano Zampini ierr = PetscBTCreate(trueNumCells, &periodicC);CHKERRQ(ierr); 511f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 512f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 513437bdd5bSStefano Zampini PetscBool pc = PETSC_FALSE; 514437bdd5bSStefano Zampini PetscInt corner; 515437bdd5bSStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 5163548e9b6SStefano Zampini pc = (PetscBool)(pc || PetscBTLookup(periodicV, gmsh_elem[c].nodes[corner] - 1)); 517437bdd5bSStefano Zampini } 518437bdd5bSStefano Zampini if (pc) { 519437bdd5bSStefano Zampini ierr = PetscBTSet(periodicC, cell);CHKERRQ(ierr); 520f45c9589SStefano Zampini ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr); 521f45c9589SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr); 5226fbe17bfSStefano Zampini } 523f45c9589SStefano Zampini cell++; 524f45c9589SStefano Zampini } 525f45c9589SStefano Zampini } 526f45c9589SStefano Zampini } 527331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 528331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 5298b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 530331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 531331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 5321d64f2ccSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, embedDim);CHKERRQ(ierr); 533331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 534331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 535331830f3SMatthew G. Knepley if (!rank) { 536f45c9589SStefano Zampini 537f45c9589SStefano Zampini if (periodicMap) { 538f45c9589SStefano Zampini PetscInt off; 539f45c9589SStefano Zampini 540f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 541f45c9589SStefano Zampini PetscInt pcone[8], corner; 542f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 5436fbe17bfSStefano Zampini if (PetscUnlikely(PetscBTLookup(periodicC,cell))) { 544f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 545f45c9589SStefano Zampini pcone[corner] = gmsh_elem[c].nodes[corner] - 1; 546f45c9589SStefano Zampini } 547f45c9589SStefano Zampini if (dim == 3) { 548f45c9589SStefano Zampini /* Tetrahedra are inverted */ 549f45c9589SStefano Zampini if (gmsh_elem[c].numNodes == 4) { 550f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 551f45c9589SStefano Zampini pcone[0] = pcone[1]; 552f45c9589SStefano Zampini pcone[1] = tmp; 553f45c9589SStefano Zampini } 554f45c9589SStefano Zampini /* Hexahedra are inverted */ 555f45c9589SStefano Zampini if (gmsh_elem[c].numNodes == 8) { 556f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 557f45c9589SStefano Zampini pcone[1] = pcone[3]; 558f45c9589SStefano Zampini pcone[3] = tmp; 559f45c9589SStefano Zampini } 560f45c9589SStefano Zampini } 561f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr); 562f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 563f45c9589SStefano Zampini PetscInt v = pcone[corner]; 564*dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 565f45c9589SStefano Zampini coords[off++] = coordsIn[v*3+d]; 566f45c9589SStefano Zampini } 567f45c9589SStefano Zampini } 5686fbe17bfSStefano Zampini } 569f45c9589SStefano Zampini cell++; 570f45c9589SStefano Zampini } 571f45c9589SStefano Zampini } 572f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 573f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 574*dd169d64SMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 575f45c9589SStefano Zampini coords[off+d] = coordsIn[periodicMapI[v]*3+d]; 576f45c9589SStefano Zampini } 577f45c9589SStefano Zampini } 578f45c9589SStefano Zampini } else { 579331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 5801d64f2ccSMatthew G. Knepley for (d = 0; d < embedDim; ++d) { 5811d64f2ccSMatthew G. Knepley coords[v*embedDim+d] = coordsIn[v*3+d]; 582331830f3SMatthew G. Knepley } 583331830f3SMatthew G. Knepley } 584331830f3SMatthew G. Knepley } 585f45c9589SStefano Zampini } 58690b157c4SStefano Zampini ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 587331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 588331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 589331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 590331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 591d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 592fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 5936fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicV);CHKERRQ(ierr); 5946fbe17bfSStefano Zampini ierr = PetscBTDestroy(&periodicC);CHKERRQ(ierr); 595a4bb7517SMichael Lange /* Clean up intermediate storage */ 59629403c87SMichael Lange if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 5973b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 598331830f3SMatthew G. Knepley PetscFunctionReturn(0); 599331830f3SMatthew G. Knepley } 600