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; 189331830f3SMatthew G. Knepley PetscScalar *coords, *coordsIn = NULL; 190f45c9589SStefano Zampini PetscInt dim = 0, coordSize, c, v, d, r, cell, *periodicMap = NULL, *periodicMapI = NULL; 191dc0ede3bSMatthew G. Knepley int i, numVertices = 0, numCells = 0, trueNumCells = 0, numRegions = 0, snum; 192331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 193331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 194*d3f73514SStefano Zampini PetscBool match, binary, bswap = PETSC_FALSE, periodic = PETSC_FALSE, usemarker = PETSC_FALSE; 195331830f3SMatthew G. Knepley PetscErrorCode ierr; 196331830f3SMatthew G. Knepley 197331830f3SMatthew G. Knepley PetscFunctionBegin; 198331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 199331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 200331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 201331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 2023b3bc66dSMichael Lange ierr = PetscLogEventBegin(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 20304d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 20404d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 205*d3f73514SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_periodic",&periodic,NULL);CHKERRQ(ierr); 206*d3f73514SStefano Zampini ierr = PetscOptionsGetBool(NULL,NULL,"-dm_plex_gmsh_usemarker",&usemarker,NULL);CHKERRQ(ierr); 207abc86ac4SMichael Lange if (!rank || binary) { 208331830f3SMatthew G. Knepley PetscBool match; 209331830f3SMatthew G. Knepley int fileType, dataSize; 210f6ac7a6aSMichael Lange float version; 211331830f3SMatthew G. Knepley 212331830f3SMatthew G. Knepley /* Read header */ 213060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 21404d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", 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"); 216060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 217f6ac7a6aSMichael Lange snum = sscanf(line, "%f %d %d", &version, &fileType, &dataSize); 218f6ac7a6aSMichael Lange if (snum != 3) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse Gmsh file header: %s", line); 219f6ac7a6aSMichael Lange if (version < 2.0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Gmsh file must be at least version 2.0"); 220331830f3SMatthew 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); 22104d1ad83SMichael Lange if (binary) { 222b9eae255SMichael Lange int checkInt; 223060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &checkInt, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 22404d1ad83SMichael Lange if (checkInt != 1) { 225b9eae255SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_ENUM, 1);CHKERRQ(ierr); 22604d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 22704d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 22804d1ad83SMichael Lange } 2290877b519SMichael Lange } else if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 230060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 23104d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 232331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 233dc0ede3bSMatthew G. Knepley /* OPTIONAL Read physical names */ 234060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 235dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$PhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 236dc0ede3bSMatthew G. Knepley if (match) { 237dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 238dc0ede3bSMatthew G. Knepley snum = sscanf(line, "%d", &numRegions); 239dc0ede3bSMatthew G. Knepley if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 240a8667449SSander Arens for (r = 0; r < numRegions; ++r) {ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr);} 241dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 242dc0ede3bSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndPhysicalNames", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 243dc0ede3bSMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 244dc0ede3bSMatthew G. Knepley /* Initial read for vertex section */ 245dc0ede3bSMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 246dc0ede3bSMatthew G. Knepley } 247dc0ede3bSMatthew G. Knepley /* Read vertices */ 24804d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 249331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 250060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 2510877b519SMichael Lange snum = sscanf(line, "%d", &numVertices); 2520877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 253854ce69bSBarry Smith ierr = PetscMalloc1(numVertices*3, &coordsIn);CHKERRQ(ierr); 25413aecfe5SMichael Lange if (binary) { 25513aecfe5SMichael Lange size_t doubleSize, intSize; 25613aecfe5SMichael Lange PetscInt elementSize; 25713aecfe5SMichael Lange char *buffer; 25813aecfe5SMichael Lange PetscScalar *baseptr; 25913aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_ENUM, &intSize);CHKERRQ(ierr); 26013aecfe5SMichael Lange ierr = PetscDataTypeGetSize(PETSC_DOUBLE, &doubleSize);CHKERRQ(ierr); 26113aecfe5SMichael Lange elementSize = (intSize + 3*doubleSize); 26213aecfe5SMichael Lange ierr = PetscMalloc1(elementSize*numVertices, &buffer);CHKERRQ(ierr); 26313aecfe5SMichael Lange ierr = PetscViewerRead(viewer, buffer, elementSize*numVertices, NULL, PETSC_CHAR);CHKERRQ(ierr); 26413aecfe5SMichael Lange if (bswap) ierr = PetscByteSwap(buffer, PETSC_CHAR, elementSize*numVertices);CHKERRQ(ierr); 265331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 26613aecfe5SMichael Lange baseptr = ((PetscScalar*)(buffer+v*elementSize+intSize)); 26713aecfe5SMichael Lange coordsIn[v*3+0] = baseptr[0]; 26813aecfe5SMichael Lange coordsIn[v*3+1] = baseptr[1]; 26913aecfe5SMichael Lange coordsIn[v*3+2] = baseptr[2]; 27013aecfe5SMichael Lange } 27113aecfe5SMichael Lange ierr = PetscFree(buffer);CHKERRQ(ierr); 27213aecfe5SMichael Lange } else { 27313aecfe5SMichael Lange for (v = 0; v < numVertices; ++v) { 274060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &i, 1, NULL, PETSC_ENUM);CHKERRQ(ierr); 275060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, NULL, PETSC_DOUBLE);CHKERRQ(ierr); 276b9eae255SMichael Lange if (i != (int)v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 277331830f3SMatthew G. Knepley } 27813aecfe5SMichael Lange } 279060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 28004d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 281331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 282331830f3SMatthew G. Knepley /* Read cells */ 283060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 28404d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", 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"); 286060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr);; 2870877b519SMichael Lange snum = sscanf(line, "%d", &numCells); 2880877b519SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 289331830f3SMatthew G. Knepley } 290db397164SMichael Lange 291abc86ac4SMichael Lange if (!rank || binary) { 292a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 293a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 294a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 295a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 2963d51f2daSMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, numCells, binary, bswap, &gmsh_elem);CHKERRQ(ierr); 297a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 298a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 299a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 300db397164SMichael Lange } 301d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 3021b82c3ebSMatthew G. Knepley ierr = PetscStrncmp(line, "$EndElements", 12, &match);CHKERRQ(ierr); 303331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 304d08df55aSStefano Zampini if (periodic) { 305f45c9589SStefano Zampini PetscInt pVert, *periodicMapT, *aux; 306d08df55aSStefano Zampini int numPeriodic; 307d08df55aSStefano Zampini 308d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 309d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$Periodic", 9, &match);CHKERRQ(ierr); 310d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 311f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapT);CHKERRQ(ierr); 312f45c9589SStefano Zampini for (i = 0; i < numVertices; i++) periodicMapT[i] = i; 313d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 314d08df55aSStefano Zampini snum = sscanf(line, "%d", &numPeriodic); 315d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 316d08df55aSStefano Zampini for (i = 0; i < numPeriodic; i++) { 317d08df55aSStefano Zampini int j, edim = -1, slaveTag = -1, masterTag = -1, nNodes; 318d08df55aSStefano Zampini 319d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 3, NULL, PETSC_STRING);CHKERRQ(ierr); 320d08df55aSStefano Zampini snum = sscanf(line, "%d %d %d", &edim, &slaveTag, &masterTag); /* slaveTag and masterTag are unused */ 321d08df55aSStefano Zampini if (snum != 3) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 322d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 323d08df55aSStefano Zampini snum = sscanf(line, "%d", &nNodes); 324d08df55aSStefano Zampini if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 325d08df55aSStefano Zampini for (j = 0; j < nNodes; j++) { 326d08df55aSStefano Zampini PetscInt slaveNode, masterNode; 327d08df55aSStefano Zampini 328d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING);CHKERRQ(ierr); 329d08df55aSStefano Zampini snum = sscanf(line, "%d %d", &slaveNode, &masterNode); 330d08df55aSStefano Zampini if (snum != 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 331f45c9589SStefano Zampini periodicMapT[slaveNode - 1] = masterNode - 1; 332d08df55aSStefano Zampini } 333d08df55aSStefano Zampini } 334d08df55aSStefano Zampini ierr = PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 335d08df55aSStefano Zampini ierr = PetscStrncmp(line, "$EndPeriodic", 12, &match);CHKERRQ(ierr); 336d08df55aSStefano Zampini if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 337d08df55aSStefano Zampini /* we may have slaves of slaves */ 338d08df55aSStefano Zampini for (i = 0; i < numVertices; i++) { 339f45c9589SStefano Zampini while (periodicMapT[periodicMapT[i]] != periodicMapT[i]) { 340f45c9589SStefano Zampini periodicMapT[i] = periodicMapT[periodicMapT[i]]; 341d08df55aSStefano Zampini } 342d08df55aSStefano Zampini } 343f45c9589SStefano Zampini /* periodicMap : from old to new numbering (periodic vertices excluded) 344f45c9589SStefano Zampini periodicMapI: from new to old numbering */ 345f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMap);CHKERRQ(ierr); 346f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &periodicMapI);CHKERRQ(ierr); 347f45c9589SStefano Zampini ierr = PetscMalloc1(numVertices, &aux);CHKERRQ(ierr); 348f45c9589SStefano Zampini for (i = 0, pVert = 0; i < numVertices; i++) { 349f45c9589SStefano Zampini if (periodicMapT[i] != i) { 350f45c9589SStefano Zampini pVert++; 351f45c9589SStefano Zampini } else { 352f45c9589SStefano Zampini aux[i] = i - pVert; 353f45c9589SStefano Zampini periodicMapI[i - pVert] = i; 354f45c9589SStefano Zampini } 355f45c9589SStefano Zampini } 356f45c9589SStefano Zampini for (i = 0 ; i < numVertices; i++) { 357f45c9589SStefano Zampini periodicMap[i] = aux[periodicMapT[i]]; 358f45c9589SStefano Zampini } 359f45c9589SStefano Zampini ierr = PetscFree(periodicMapT);CHKERRQ(ierr); 360f45c9589SStefano Zampini ierr = PetscFree(aux);CHKERRQ(ierr); 361f45c9589SStefano Zampini /* remove periodic vertices */ 362f45c9589SStefano Zampini numVertices = numVertices - pVert; 363d08df55aSStefano Zampini } 364331830f3SMatthew G. Knepley } 365abc86ac4SMichael Lange /* For binary we read on all ranks, but only build the plex on rank 0 */ 366abc86ac4SMichael Lange if (binary && rank) {trueNumCells = 0; numVertices = 0;}; 367a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 368331830f3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 369331830f3SMatthew G. Knepley if (!rank) { 370a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 371a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 372a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 373a4bb7517SMichael Lange cell++; 374331830f3SMatthew G. Knepley } 375331830f3SMatthew G. Knepley } 376331830f3SMatthew G. Knepley } 377331830f3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 378a4bb7517SMichael Lange /* Add cell-vertex connections */ 379331830f3SMatthew G. Knepley if (!rank) { 380331830f3SMatthew G. Knepley PetscInt pcone[8], corner; 381a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 382a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 383a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 384d08df55aSStefano Zampini pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + trueNumCells 385d08df55aSStefano Zampini : gmsh_elem[c].nodes[corner] -1 + trueNumCells; 386331830f3SMatthew G. Knepley } 38797ed6be6Ssarens if (dim == 3) { 38897ed6be6Ssarens /* Tetrahedra are inverted */ 38997ed6be6Ssarens if (gmsh_elem[c].numNodes == 4) { 39097ed6be6Ssarens PetscInt tmp = pcone[0]; 39197ed6be6Ssarens pcone[0] = pcone[1]; 39297ed6be6Ssarens pcone[1] = tmp; 39397ed6be6Ssarens } 39497ed6be6Ssarens /* Hexahedra are inverted */ 39597ed6be6Ssarens if (gmsh_elem[c].numNodes == 8) { 39697ed6be6Ssarens PetscInt tmp = pcone[1]; 39797ed6be6Ssarens pcone[1] = pcone[3]; 39897ed6be6Ssarens pcone[3] = tmp; 39997ed6be6Ssarens } 40097ed6be6Ssarens } 401a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 402a4bb7517SMichael Lange cell++; 403331830f3SMatthew G. Knepley } 404a4bb7517SMichael Lange } 405331830f3SMatthew G. Knepley } 4066228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 407c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 408331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 409331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 410331830f3SMatthew G. Knepley if (interpolate) { 411e56d480eSMatthew G. Knepley DM idm = NULL; 412331830f3SMatthew G. Knepley 413331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 414331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 415331830f3SMatthew G. Knepley *dm = idm; 416331830f3SMatthew G. Knepley } 417*d3f73514SStefano Zampini if (usemarker) { 418*d3f73514SStefano Zampini PetscInt f, fStart, fEnd; 419*d3f73514SStefano Zampini 420*d3f73514SStefano Zampini ierr = DMCreateLabel(*dm,"marker");CHKERRQ(ierr); 421*d3f73514SStefano Zampini ierr = DMPlexGetHeightStratum(*dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 422*d3f73514SStefano Zampini for (f = fStart; f < fEnd; ++f) { 423*d3f73514SStefano Zampini PetscInt suppSize; 424*d3f73514SStefano Zampini 425*d3f73514SStefano Zampini ierr = DMPlexGetSupportSize(*dm, f, &suppSize);CHKERRQ(ierr); 426*d3f73514SStefano Zampini if (suppSize == 1) { 427*d3f73514SStefano Zampini PetscInt *cone = NULL, coneSize, p; 428*d3f73514SStefano Zampini 429*d3f73514SStefano Zampini ierr = DMPlexGetTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 430*d3f73514SStefano Zampini for (p = 0; p < coneSize; p += 2) { 431*d3f73514SStefano Zampini ierr = DMSetLabelValue(*dm, "marker", cone[p], 1);CHKERRQ(ierr); 432*d3f73514SStefano Zampini } 433*d3f73514SStefano Zampini ierr = DMPlexRestoreTransitiveClosure(*dm, f, PETSC_TRUE, &coneSize, &cone);CHKERRQ(ierr); 434*d3f73514SStefano Zampini } 435*d3f73514SStefano Zampini } 436*d3f73514SStefano Zampini } 43716ad7e67SMichael Lange 43816ad7e67SMichael Lange if (!rank) { 43916ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 44016ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 441d1a54cefSMatthew G. Knepley 44216ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 44316ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 444a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 44516ad7e67SMichael Lange PetscInt joinSize; 44616ad7e67SMichael Lange const PetscInt *join; 447a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 448d08df55aSStefano Zampini pcone[corner] = periodicMap ? periodicMap[gmsh_elem[c].nodes[corner] - 1] + vStart 449d08df55aSStefano Zampini : gmsh_elem[c].nodes[corner] - 1 + vStart; 45016ad7e67SMichael Lange } 451a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 452a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 453c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 454a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 45516ad7e67SMichael Lange } 456a4bb7517SMichael Lange } 4576e1dd89aSLawrence Mitchell 4586e1dd89aSLawrence Mitchell /* Create cell sets */ 4596e1dd89aSLawrence Mitchell for (cell = 0, c = 0; c < numCells; ++c) { 4606e1dd89aSLawrence Mitchell if (gmsh_elem[c].dim == dim) { 4616e1dd89aSLawrence Mitchell if (gmsh_elem[c].numTags > 0) { 4626e1dd89aSLawrence Mitchell ierr = DMSetLabelValue(*dm, "Cell Sets", cell, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 4636e1dd89aSLawrence Mitchell cell++; 4646e1dd89aSLawrence Mitchell } 4656e1dd89aSLawrence Mitchell } 4666e1dd89aSLawrence Mitchell } 4670c070f12SSander Arens 4680c070f12SSander Arens /* Create vertex sets */ 4690c070f12SSander Arens for (c = 0; c < numCells; ++c) { 4700c070f12SSander Arens if (gmsh_elem[c].dim == 0) { 4710c070f12SSander Arens if (gmsh_elem[c].numTags > 0) { 472d08df55aSStefano Zampini PetscInt vid = periodicMap ? periodicMap[gmsh_elem[c].nodes[0] - 1] + vStart 473d08df55aSStefano Zampini : gmsh_elem[c].nodes[0] - 1 + vStart; 474d08df55aSStefano Zampini ierr = DMSetLabelValue(*dm, "Vertex Sets", vid, gmsh_elem[c].tags[0]);CHKERRQ(ierr); 4750c070f12SSander Arens } 4760c070f12SSander Arens } 4770c070f12SSander Arens } 47816ad7e67SMichael Lange } 47916ad7e67SMichael Lange 480331830f3SMatthew G. Knepley /* Read coordinates */ 481979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 482331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 483331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 484f45c9589SStefano Zampini if (periodicMap) { /* we need to localize coordinates on cells */ 485f45c9589SStefano Zampini ierr = PetscSectionSetChart(coordSection, 0, trueNumCells + numVertices);CHKERRQ(ierr); 486f45c9589SStefano Zampini } else { 48775b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 488f45c9589SStefano Zampini } 48975b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 490331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 491331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 492331830f3SMatthew G. Knepley } 493f45c9589SStefano Zampini if (periodicMap) { 494f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 495f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 496f45c9589SStefano Zampini ierr = PetscSectionSetDof(coordSection, cell, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr); 497f45c9589SStefano Zampini ierr = PetscSectionSetFieldDof(coordSection, cell, 0, gmsh_elem[c].numNodes*dim);CHKERRQ(ierr); 498f45c9589SStefano Zampini cell++; 499f45c9589SStefano Zampini } 500f45c9589SStefano Zampini } 501f45c9589SStefano Zampini } 502331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 503331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 5048b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 505331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 506331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 5077635fff5Ssarens ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 508331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 509331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 510331830f3SMatthew G. Knepley if (!rank) { 511f45c9589SStefano Zampini 512f45c9589SStefano Zampini if (periodicMap) { 513f45c9589SStefano Zampini PetscInt off; 514f45c9589SStefano Zampini 515f45c9589SStefano Zampini for (cell = 0, c = 0; c < numCells; ++c) { 516f45c9589SStefano Zampini PetscInt pcone[8], corner; 517f45c9589SStefano Zampini if (gmsh_elem[c].dim == dim) { 518f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 519f45c9589SStefano Zampini pcone[corner] = gmsh_elem[c].nodes[corner] - 1; 520f45c9589SStefano Zampini } 521f45c9589SStefano Zampini if (dim == 3) { 522f45c9589SStefano Zampini /* Tetrahedra are inverted */ 523f45c9589SStefano Zampini if (gmsh_elem[c].numNodes == 4) { 524f45c9589SStefano Zampini PetscInt tmp = pcone[0]; 525f45c9589SStefano Zampini pcone[0] = pcone[1]; 526f45c9589SStefano Zampini pcone[1] = tmp; 527f45c9589SStefano Zampini } 528f45c9589SStefano Zampini /* Hexahedra are inverted */ 529f45c9589SStefano Zampini if (gmsh_elem[c].numNodes == 8) { 530f45c9589SStefano Zampini PetscInt tmp = pcone[1]; 531f45c9589SStefano Zampini pcone[1] = pcone[3]; 532f45c9589SStefano Zampini pcone[3] = tmp; 533f45c9589SStefano Zampini } 534f45c9589SStefano Zampini } 535f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, cell, &off);CHKERRQ(ierr); 536f45c9589SStefano Zampini for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 537f45c9589SStefano Zampini PetscInt v = pcone[corner]; 538f45c9589SStefano Zampini for (d = 0; d < dim; ++d) { 539f45c9589SStefano Zampini coords[off++] = coordsIn[v*3+d]; 540f45c9589SStefano Zampini } 541f45c9589SStefano Zampini } 542f45c9589SStefano Zampini cell++; 543f45c9589SStefano Zampini } 544f45c9589SStefano Zampini } 545f45c9589SStefano Zampini for (v = 0; v < numVertices; ++v) { 546f45c9589SStefano Zampini ierr = PetscSectionGetOffset(coordSection, v + trueNumCells, &off);CHKERRQ(ierr); 547f45c9589SStefano Zampini for (d = 0; d < dim; ++d) { 548f45c9589SStefano Zampini coords[off+d] = coordsIn[periodicMapI[v]*3+d]; 549f45c9589SStefano Zampini } 550f45c9589SStefano Zampini } 551f45c9589SStefano Zampini } else { 552331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 553331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 554331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 555331830f3SMatthew G. Knepley } 556331830f3SMatthew G. Knepley } 557331830f3SMatthew G. Knepley } 558f45c9589SStefano Zampini } 55990b157c4SStefano Zampini ierr = DMSetPeriodicity(*dm, periodic, NULL, NULL, NULL);CHKERRQ(ierr); 560331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 561331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 562331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 563331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 564d08df55aSStefano Zampini ierr = PetscFree(periodicMap);CHKERRQ(ierr); 565fcd9ca0aSStefano Zampini ierr = PetscFree(periodicMapI);CHKERRQ(ierr); 566a4bb7517SMichael Lange /* Clean up intermediate storage */ 56729403c87SMichael Lange if (!rank || binary) ierr = PetscFree(gmsh_elem);CHKERRQ(ierr); 5683b3bc66dSMichael Lange ierr = PetscLogEventEnd(DMPLEX_CreateGmsh,*dm,0,0,0);CHKERRQ(ierr); 569331830f3SMatthew G. Knepley PetscFunctionReturn(0); 570331830f3SMatthew G. Knepley } 571