1f9720913SMatthew G. Knepley #define PETSCDM_DLL 2f9720913SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3f9720913SMatthew G. Knepley 4f9720913SMatthew G. Knepley /*@C 5f9720913SMatthew G. Knepley DMPlexCreatePLYFromFile - Create a DMPlex mesh from a PLY file. 6f9720913SMatthew G. Knepley 7f9720913SMatthew G. Knepley + comm - The MPI communicator 8f9720913SMatthew G. Knepley . filename - Name of the .med file 9f9720913SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 10f9720913SMatthew G. Knepley 11f9720913SMatthew G. Knepley Output Parameter: 12f9720913SMatthew G. Knepley . dm - The DM object representing the mesh 13f9720913SMatthew G. Knepley 14f9720913SMatthew G. Knepley Note: https://en.wikipedia.org/wiki/PLY_(file_format) 15f9720913SMatthew G. Knepley 16f9720913SMatthew G. Knepley Level: beginner 17f9720913SMatthew G. Knepley 18f9720913SMatthew G. Knepley .seealso: DMPlexCreateFromFile(), DMPlexCreateMedFromFile(), DMPlexCreateGmsh(), DMPlexCreate() 19f9720913SMatthew G. Knepley @*/ 20f9720913SMatthew G. Knepley PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 21f9720913SMatthew G. Knepley { 22f9720913SMatthew G. Knepley PetscViewer viewer; 23f9720913SMatthew G. Knepley Vec coordinates; 24f9720913SMatthew G. Knepley PetscSection coordSection; 25f9720913SMatthew G. Knepley PetscScalar *coords; 26f9720913SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN], ntype[16], itype[16], name[1024], vtype[16]; 27f9720913SMatthew G. Knepley PetscBool match, byteSwap = PETSC_FALSE; 2864b21ff4SMatthew G. Knepley PetscInt dim = 2, cdim = 3, Nvp = 0, coordSize, xi = -1, yi = -1, zi = -1, v, c, p; 29f9720913SMatthew G. Knepley PetscMPIInt rank; 30f9720913SMatthew G. Knepley int snum, Nv, Nc; 31f9720913SMatthew G. Knepley 32f9720913SMatthew G. Knepley PetscFunctionBegin; 335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(comm, &viewer)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer, PETSCVIEWERBINARY)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer, filename)); 38dd400576SPatrick Sanan if (rank == 0) { 39f9720913SMatthew G. Knepley PetscBool isAscii, isBinaryBig, isBinaryLittle; 40f9720913SMatthew G. Knepley 41f9720913SMatthew G. Knepley /* Check for PLY file */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match)); 44*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 45f9720913SMatthew G. Knepley /* Check PLY format */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match)); 48*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle)); 53*28b400f6SJacob Faibussowitsch PetscCheck(!isAscii,PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported"); 54f9720913SMatthew G. Knepley else if (isBinaryLittle) byteSwap = PETSC_TRUE; 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match)); 57*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid version of PLY file, %s", line); 58f9720913SMatthew G. Knepley /* Ignore comments */ 59f9720913SMatthew G. Knepley match = PETSC_TRUE; 60f9720913SMatthew G. Knepley while (match) { 61f9720913SMatthew G. Knepley char buf = '\0'; 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match)); 645f80ce2aSJacob Faibussowitsch if (match) while (buf != '\n') CHKERRQ(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR)); 65f9720913SMatthew G. Knepley } 66f9720913SMatthew G. Knepley /* Read vertex information */ 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 68*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match)); 71*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 73f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nv); 742c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 75f9720913SMatthew G. Knepley match = PETSC_TRUE; 76f9720913SMatthew G. Knepley while (match) { 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match)); 79f9720913SMatthew G. Knepley if (match) { 80f9720913SMatthew G. Knepley PetscBool matchB; 81f9720913SMatthew G. Knepley 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 832c71b3e2SJacob Faibussowitsch PetscCheckFalse(Nvp >= 16,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle more than 16 property statements in PLY file header: %s", line); 84f9720913SMatthew G. Knepley snum = sscanf(line, "%s %s", ntype, name); 852c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(ntype, "float32", 16, &matchB)); 87f9720913SMatthew G. Knepley if (matchB) { 88f9720913SMatthew G. Knepley vtype[Nvp] = 'f'; 89f9720913SMatthew G. Knepley } else { 905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(ntype, "int32", 16, &matchB)); 91f9720913SMatthew G. Knepley if (matchB) { 92f9720913SMatthew G. Knepley vtype[Nvp] = 'd'; 93f9720913SMatthew G. Knepley } else { 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(ntype, "uint8", 16, &matchB)); 95f9720913SMatthew G. Knepley if (matchB) { 96f9720913SMatthew G. Knepley vtype[Nvp] = 'c'; 9798921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse type in PLY file header: %s", line); 98f9720913SMatthew G. Knepley } 99f9720913SMatthew G. Knepley } 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(name, "x", 16, &matchB)); 101f9720913SMatthew G. Knepley if (matchB) {xi = Nvp;} 1025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(name, "y", 16, &matchB)); 103f9720913SMatthew G. Knepley if (matchB) {yi = Nvp;} 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(name, "z", 16, &matchB)); 105f9720913SMatthew G. Knepley if (matchB) {zi = Nvp;} 106f9720913SMatthew G. Knepley ++Nvp; 107f9720913SMatthew G. Knepley } 108f9720913SMatthew G. Knepley } 109f9720913SMatthew G. Knepley /* Read cell information */ 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 111*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "face", PETSC_MAX_PATH_LEN, &match)); 114*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 116f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nc); 1172c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING)); 119f9720913SMatthew G. Knepley snum = sscanf(line, "property list %s %s %s", ntype, itype, name); 1202c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 3,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(ntype, "uint8", 1024, &match)); 122*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(name, "vertex_indices", 1024, &match)); 124*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line); 125f9720913SMatthew G. Knepley /* Header should terminate */ 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match)); 128*28b400f6SJacob Faibussowitsch PetscCheck(match,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file: %s", line); 12925ce1634SJed Brown } else { 13025ce1634SJed Brown Nc = Nv = 0; 131f9720913SMatthew G. Knepley } 1325f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 1335f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetChart(*dm, 0, Nc+Nv)); 1355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinateDim(*dm, cdim)); 137f9720913SMatthew G. Knepley /* Read coordinates */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(*dm, &coordSection)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(coordSection, 1)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(coordSection, Nc, Nc + Nv)); 142f9720913SMatthew G. Knepley for (v = Nc; v < Nc+Nv; ++v) { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(coordSection, v, cdim)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 145f9720913SMatthew G. Knepley } 1465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(coordSection)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetStorageSize(coordSection, &coordSize)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF, &coordinates)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(coordinates, cdim)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetType(coordinates, VECSTANDARD)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 154dd400576SPatrick Sanan if (rank == 0) { 155f9720913SMatthew G. Knepley float rbuf[1]; 156f9720913SMatthew G. Knepley int ibuf[1]; 157f9720913SMatthew G. Knepley 158f9720913SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 159f9720913SMatthew G. Knepley for (p = 0; p < Nvp; ++p) { 160f9720913SMatthew G. Knepley if (vtype[p] == 'f') { 1615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT)); 1625f80ce2aSJacob Faibussowitsch if (byteSwap) CHKERRQ(PetscByteSwap(&rbuf, PETSC_FLOAT, 1)); 163f9720913SMatthew G. Knepley if (p == xi) coords[v*cdim+0] = rbuf[0]; 164f9720913SMatthew G. Knepley else if (p == yi) coords[v*cdim+1] = rbuf[0]; 165f9720913SMatthew G. Knepley else if (p == zi) coords[v*cdim+2] = rbuf[0]; 166f9720913SMatthew G. Knepley } else if (vtype[p] == 'd') { 1675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT)); 1685f80ce2aSJacob Faibussowitsch if (byteSwap) CHKERRQ(PetscByteSwap(&ibuf, PETSC_INT, 1)); 169f9720913SMatthew G. Knepley } else if (vtype[p] == 'c') { 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 171f9720913SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid vertex property type in PLY file"); 172f9720913SMatthew G. Knepley } 173f9720913SMatthew G. Knepley } 174f9720913SMatthew G. Knepley } 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinatesLocal(*dm, coordinates)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&coordinates)); 178f9720913SMatthew G. Knepley /* Read topology */ 179dd400576SPatrick Sanan if (rank == 0) { 180f9720913SMatthew G. Knepley char ibuf[1]; 181f9720913SMatthew G. Knepley PetscInt vbuf[16], corners; 182f9720913SMatthew G. Knepley 183f9720913SMatthew G. Knepley /* Assume same cells */ 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 185f9720913SMatthew G. Knepley corners = ibuf[0]; 1865f80ce2aSJacob Faibussowitsch for (c = 0; c < Nc; ++c) CHKERRQ(DMPlexSetConeSize(*dm, c, corners)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(*dm)); 188f9720913SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 189f9720913SMatthew G. Knepley if (c > 0) { 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 191f9720913SMatthew G. Knepley } 1922c71b3e2SJacob Faibussowitsch PetscCheckFalse(ibuf[0] != corners,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "All cells must have the same number of vertices in PLY file: %D != %D", ibuf[0], corners); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT)); 1945f80ce2aSJacob Faibussowitsch if (byteSwap) CHKERRQ(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0])); 195f9720913SMatthew G. Knepley for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc; 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCone(*dm, c, vbuf)); 197f9720913SMatthew G. Knepley } 198f9720913SMatthew G. Knepley } 1995f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSymmetrize(*dm)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexStratify(*dm)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 202f9720913SMatthew G. Knepley if (interpolate) { 2035fd9971aSMatthew G. Knepley DM idm; 204f9720913SMatthew G. Knepley 2055f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 207f9720913SMatthew G. Knepley *dm = idm; 208f9720913SMatthew G. Knepley } 209f9720913SMatthew G. Knepley PetscFunctionReturn(0); 210f9720913SMatthew G. Knepley } 211