1331830f3SMatthew G. Knepley #define PETSCDM_DLL 2331830f3SMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3331830f3SMatthew G. Knepley 4331830f3SMatthew G. Knepley #undef __FUNCT__ 5331830f3SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateGmsh" 6331830f3SMatthew G. Knepley /*@ 7331830f3SMatthew G. Knepley DMPlexCreateGmsh - Create a DMPlex mesh from a Gmsh file. 8331830f3SMatthew G. Knepley 9331830f3SMatthew G. Knepley Collective on comm 10331830f3SMatthew G. Knepley 11331830f3SMatthew G. Knepley Input Parameters: 12331830f3SMatthew G. Knepley + comm - The MPI communicator 13331830f3SMatthew G. Knepley . viewer - The Viewer associated with a Gmsh file 14331830f3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 15331830f3SMatthew G. Knepley 16331830f3SMatthew G. Knepley Output Parameter: 17331830f3SMatthew G. Knepley . dm - The DM object representing the mesh 18331830f3SMatthew G. Knepley 19331830f3SMatthew G. Knepley Note: http://www.geuz.org/gmsh/doc/texinfo/#MSH-ASCII-file-format 20331830f3SMatthew G. Knepley 21331830f3SMatthew G. Knepley Level: beginner 22331830f3SMatthew G. Knepley 23331830f3SMatthew G. Knepley .keywords: mesh,Gmsh 24331830f3SMatthew G. Knepley .seealso: DMPLEX, DMCreate() 25331830f3SMatthew G. Knepley @*/ 26331830f3SMatthew G. Knepley PetscErrorCode DMPlexCreateGmsh(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 27331830f3SMatthew G. Knepley { 28*04d1ad83SMichael Lange PetscViewerType vtype; 29a4bb7517SMichael Lange GmshElement *gmsh_elem; 30331830f3SMatthew G. Knepley PetscSection coordSection; 31331830f3SMatthew G. Knepley Vec coordinates; 32331830f3SMatthew G. Knepley PetscScalar *coords, *coordsIn = NULL; 33a4bb7517SMichael Lange PetscInt dim = 0, coordSize, c, v, d, cell; 34a4bb7517SMichael Lange int numVertices = 0, numCells = 0, trueNumCells = 0, snum; 35331830f3SMatthew G. Knepley PetscMPIInt num_proc, rank; 36331830f3SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN]; 37*04d1ad83SMichael Lange PetscBool match, binary, bswap = PETSC_FALSE; 38331830f3SMatthew G. Knepley PetscErrorCode ierr; 39331830f3SMatthew G. Knepley 40331830f3SMatthew G. Knepley PetscFunctionBegin; 41331830f3SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 42331830f3SMatthew G. Knepley ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 43331830f3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 44331830f3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 45*04d1ad83SMichael Lange ierr = PetscViewerGetType(viewer, &vtype);CHKERRQ(ierr); 46*04d1ad83SMichael Lange ierr = PetscStrcmp(vtype, PETSCVIEWERBINARY, &binary);CHKERRQ(ierr); 47331830f3SMatthew G. Knepley if (!rank) { 48331830f3SMatthew G. Knepley PetscBool match; 49331830f3SMatthew G. Knepley int fileType, dataSize; 50331830f3SMatthew G. Knepley 51331830f3SMatthew G. Knepley /* Read header */ 52*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr); 53*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$MeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 54331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 55*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 3, PETSC_STRING);CHKERRQ(ierr); 56*04d1ad83SMichael Lange snum = sscanf(line, "2.2 %d %d", &fileType, &dataSize);CHKERRQ(snum != 2); 57331830f3SMatthew 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); 58*04d1ad83SMichael Lange if (binary) { 59*04d1ad83SMichael Lange PetscInt checkInt; 60*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &checkInt, 1, PETSC_INT);CHKERRQ(ierr); 61*04d1ad83SMichael Lange if (checkInt != 1) { 62*04d1ad83SMichael Lange ierr = PetscByteSwap(&checkInt, PETSC_INT, 1);CHKERRQ(ierr); 63*04d1ad83SMichael Lange if (checkInt == 1) bswap = PETSC_TRUE; 64*04d1ad83SMichael Lange else SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh binary file", fileType); 65*04d1ad83SMichael Lange } 66*04d1ad83SMichael Lange } else { 67*04d1ad83SMichael Lange if (fileType) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File type %d is not a valid Gmsh ASCII file", fileType); 68*04d1ad83SMichael Lange } 69*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr); 70*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndMeshFormat", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 71331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 72331830f3SMatthew G. Knepley /* Read vertices */ 73*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr); 74*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Nodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 75331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 76*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr); 77*04d1ad83SMichael Lange snum = sscanf(line, "%d", &numVertices);CHKERRQ(snum != 1); 78331830f3SMatthew G. Knepley ierr = PetscMalloc(numVertices*3 * sizeof(PetscScalar), &coordsIn);CHKERRQ(ierr); 79331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 80331830f3SMatthew G. Knepley int i; 81*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &i, 1, PETSC_INT);CHKERRQ(ierr); 82*04d1ad83SMichael Lange if (bswap) ierr = PetscByteSwap(&i, PETSC_INT, 1);CHKERRQ(ierr); 83*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &(coordsIn[v*3]), 3, PETSC_DOUBLE);CHKERRQ(ierr); 84*04d1ad83SMichael Lange if (bswap) ierr = PetscByteSwap(&(coordsIn[v*3]), PETSC_DOUBLE, 3);CHKERRQ(ierr); 85331830f3SMatthew G. Knepley if (i != v+1) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid node number %d should be %d", i, v+1); 86331830f3SMatthew G. Knepley } 87*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);; 88*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndNodes", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 89331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 90331830f3SMatthew G. Knepley /* Read cells */ 91*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);; 92*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$Elements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 93331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 94*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);; 95*04d1ad83SMichael Lange snum = sscanf(line, "%d", &numCells);CHKERRQ(snum != 1); 96331830f3SMatthew G. Knepley } 97db397164SMichael Lange 98331830f3SMatthew G. Knepley if (!rank) { 99a4bb7517SMichael Lange /* Gmsh elements can be of any dimension/co-dimension, so we need to traverse the 100a4bb7517SMichael Lange file contents multiple times to figure out the true number of cells and facets 101a4bb7517SMichael Lange in the given mesh. To make this more efficient we read the file contents only 102a4bb7517SMichael Lange once and store them in memory, while determining the true number of cells. */ 103a4bb7517SMichael Lange ierr = PetscMalloc1(numCells, &gmsh_elem);CHKERRQ(ierr); 104a4bb7517SMichael Lange for (trueNumCells=0, c = 0; c < numCells; ++c) { 105*04d1ad83SMichael Lange ierr = DMPlexCreateGmsh_ReadElement(viewer, binary, bswap, &gmsh_elem[c]);CHKERRQ(ierr); 106a4bb7517SMichael Lange if (gmsh_elem[c].dim > dim) {dim = gmsh_elem[c].dim; trueNumCells = 0;} 107a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) trueNumCells++; 108331830f3SMatthew G. Knepley } 109*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, line, 1, PETSC_STRING);CHKERRQ(ierr);; 110*04d1ad83SMichael Lange ierr = PetscStrncmp(line, "$EndElements", PETSC_MAX_PATH_LEN, &match);CHKERRQ(ierr); 111331830f3SMatthew G. Knepley if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Gmsh file"); 112331830f3SMatthew G. Knepley } 113a4bb7517SMichael Lange /* Allocate the cell-vertex mesh */ 114a4bb7517SMichael Lange ierr = DMPlexSetChart(*dm, 0, trueNumCells+numVertices);CHKERRQ(ierr); 115a4bb7517SMichael Lange if (!rank) { 116a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 117a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 118a4bb7517SMichael Lange ierr = DMPlexSetConeSize(*dm, cell, gmsh_elem[c].numNodes);CHKERRQ(ierr); 119a4bb7517SMichael Lange cell++; 120a4bb7517SMichael Lange } 121a4bb7517SMichael Lange } 122a4bb7517SMichael Lange } 123a4bb7517SMichael Lange ierr = DMSetUp(*dm);CHKERRQ(ierr); 124a4bb7517SMichael Lange /* Add cell-vertex connections */ 125a4bb7517SMichael Lange if (!rank) { 126a4bb7517SMichael Lange PetscInt pcone[8], corner; 127a4bb7517SMichael Lange for (cell = 0, c = 0; c < numCells; ++c) { 128a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim) { 129a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 130a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + trueNumCells-1; 131a4bb7517SMichael Lange } 132a4bb7517SMichael Lange ierr = DMPlexSetCone(*dm, cell, pcone);CHKERRQ(ierr); 133a4bb7517SMichael Lange cell++; 134a4bb7517SMichael Lange } 135a4bb7517SMichael Lange } 136a4bb7517SMichael Lange } 1376228fc9fSJed Brown ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 138a4bb7517SMichael Lange ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 139331830f3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 140331830f3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 141331830f3SMatthew G. Knepley if (interpolate) { 142e56d480eSMatthew G. Knepley DM idm = NULL; 143331830f3SMatthew G. Knepley 144331830f3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 145331830f3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 146331830f3SMatthew G. Knepley *dm = idm; 147331830f3SMatthew G. Knepley } 14816ad7e67SMichael Lange 14916ad7e67SMichael Lange if (!rank) { 15016ad7e67SMichael Lange /* Apply boundary IDs by finding the relevant facets with vertex joins */ 15116ad7e67SMichael Lange PetscInt pcone[8], corner, vStart, vEnd; 152d1a54cefSMatthew G. Knepley 15316ad7e67SMichael Lange ierr = DMPlexGetDepthStratum(*dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 15416ad7e67SMichael Lange for (c = 0; c < numCells; ++c) { 155a4bb7517SMichael Lange if (gmsh_elem[c].dim == dim-1) { 15616ad7e67SMichael Lange PetscInt joinSize; 15716ad7e67SMichael Lange const PetscInt *join; 158a4bb7517SMichael Lange for (corner = 0; corner < gmsh_elem[c].numNodes; ++corner) { 159a4bb7517SMichael Lange pcone[corner] = gmsh_elem[c].nodes[corner] + vStart - 1; 16016ad7e67SMichael Lange } 161a4bb7517SMichael Lange ierr = DMPlexGetFullJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 162a4bb7517SMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for element %d", gmsh_elem[c].id); 163a4bb7517SMichael Lange ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], gmsh_elem[c].tags[0]);CHKERRQ(ierr); 164a4bb7517SMichael Lange ierr = DMPlexRestoreJoin(*dm, gmsh_elem[c].numNodes, (const PetscInt *) pcone, &joinSize, &join);CHKERRQ(ierr); 16516ad7e67SMichael Lange } 166a4bb7517SMichael Lange } 16716ad7e67SMichael Lange } 16816ad7e67SMichael Lange 169331830f3SMatthew G. Knepley /* Read coordinates */ 170979c4b60SMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 171331830f3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 172331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 17375b5763bSMichael Lange ierr = PetscSectionSetChart(coordSection, trueNumCells, trueNumCells + numVertices);CHKERRQ(ierr); 17475b5763bSMichael Lange for (v = trueNumCells; v < trueNumCells+numVertices; ++v) { 175331830f3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 176331830f3SMatthew G. Knepley ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 177331830f3SMatthew G. Knepley } 178331830f3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 179331830f3SMatthew G. Knepley ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 180331830f3SMatthew G. Knepley ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 181331830f3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 182331830f3SMatthew G. Knepley ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 183331830f3SMatthew G. Knepley ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 184331830f3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 185331830f3SMatthew G. Knepley if (!rank) { 186331830f3SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 187331830f3SMatthew G. Knepley for (d = 0; d < dim; ++d) { 188331830f3SMatthew G. Knepley coords[v*dim+d] = coordsIn[v*3+d]; 189331830f3SMatthew G. Knepley } 190331830f3SMatthew G. Knepley } 191331830f3SMatthew G. Knepley } 192331830f3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 193331830f3SMatthew G. Knepley ierr = PetscFree(coordsIn);CHKERRQ(ierr); 194331830f3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 195331830f3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 196a4bb7517SMichael Lange /* Clean up intermediate storage */ 197a4bb7517SMichael Lange if (!rank) { 198a4bb7517SMichael Lange for (c = 0; c < numCells; ++c) { 199a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem[c].nodes); 200a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem[c].tags); 201a4bb7517SMichael Lange } 202a4bb7517SMichael Lange ierr = PetscFree(gmsh_elem); 203a4bb7517SMichael Lange } 204331830f3SMatthew G. Knepley PetscFunctionReturn(0); 205331830f3SMatthew G. Knepley } 206db397164SMichael Lange 207db397164SMichael Lange #undef __FUNCT__ 20830412aabSMichael Lange #define __FUNCT__ "DMPlexCreateGmsh_ReadElement" 209*04d1ad83SMichael Lange PetscErrorCode DMPlexCreateGmsh_ReadElement(PetscViewer viewer, PetscBool binary, PetscBool byteSwap, GmshElement *ele) 210db397164SMichael Lange { 211*04d1ad83SMichael Lange int cellType, numElem; 212a4bb7517SMichael Lange PetscErrorCode ierr; 213db397164SMichael Lange 21430412aabSMichael Lange PetscFunctionBegin; 215*04d1ad83SMichael Lange if (binary) { 216*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_INT);CHKERRQ(ierr); 217*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap(&cellType, PETSC_INT, 1);CHKERRQ(ierr); 218*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &numElem, 1, PETSC_INT);CHKERRQ(ierr); 219*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap(&numElem, PETSC_INT, 1);CHKERRQ(ierr); 220*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_INT);CHKERRQ(ierr); 221*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap(&(ele->numTags), PETSC_INT, 1);CHKERRQ(ierr); 222*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_INT);CHKERRQ(ierr); 223*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap( &(ele->id), PETSC_INT, 1);CHKERRQ(ierr); 224*04d1ad83SMichael Lange } else { 225*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &(ele->id), 1, PETSC_INT);CHKERRQ(ierr); 226*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &cellType, 1, PETSC_INT);CHKERRQ(ierr); 227*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, &(ele->numTags), 1, PETSC_INT);CHKERRQ(ierr); 228*04d1ad83SMichael Lange } 229a4bb7517SMichael Lange ierr = PetscMalloc1(ele->numTags, &(ele->tags));CHKERRQ(ierr); 230*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, ele->tags, ele->numTags, PETSC_INT);CHKERRQ(ierr); 231*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap(ele->tags, PETSC_INT, ele->numTags);CHKERRQ(ierr); 23230412aabSMichael Lange switch (cellType) { 23330412aabSMichael Lange case 1: /* 2-node line */ 234a4bb7517SMichael Lange ele->dim = 1; 235a4bb7517SMichael Lange ele->numNodes = 2; 23630412aabSMichael Lange break; 23730412aabSMichael Lange case 2: /* 3-node triangle */ 238a4bb7517SMichael Lange ele->dim = 2; 239a4bb7517SMichael Lange ele->numNodes = 3; 24030412aabSMichael Lange break; 24130412aabSMichael Lange case 3: /* 4-node quadrangle */ 242a4bb7517SMichael Lange ele->dim = 2; 243a4bb7517SMichael Lange ele->numNodes = 4; 24430412aabSMichael Lange break; 24530412aabSMichael Lange case 4: /* 4-node tetrahedron */ 246a4bb7517SMichael Lange ele->dim = 3; 247a4bb7517SMichael Lange ele->numNodes = 4; 24830412aabSMichael Lange break; 24930412aabSMichael Lange case 5: /* 8-node hexahedron */ 250a4bb7517SMichael Lange ele->dim = 3; 251a4bb7517SMichael Lange ele->numNodes = 8; 25230412aabSMichael Lange break; 2536b7f3382SJustin Chang case 15: /* 1-node vertex */ 254a4bb7517SMichael Lange ele->dim = 0; 255a4bb7517SMichael Lange ele->numNodes = 1; 2566b7f3382SJustin Chang break; 25730412aabSMichael Lange default: 25830412aabSMichael Lange SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported Gmsh element type %d", cellType); 25930412aabSMichael Lange } 260*04d1ad83SMichael Lange ierr = PetscMalloc1(ele->numNodes, &(ele->nodes));CHKERRQ(ierr); 261*04d1ad83SMichael Lange ierr = PetscViewerRead(viewer, ele->nodes, ele->numNodes, PETSC_INT);CHKERRQ(ierr); 262*04d1ad83SMichael Lange if (byteSwap) ierr = PetscByteSwap(ele->nodes, PETSC_INT, ele->numNodes);CHKERRQ(ierr); 26330412aabSMichael Lange PetscFunctionReturn(0); 264db397164SMichael Lange } 265