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 18db781477SPatrick Sanan .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 19f9720913SMatthew G. Knepley @*/ 20d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 21d71ae5a4SJacob Faibussowitsch { 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; 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 359566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 379566063dSJacob Faibussowitsch PetscCall(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 */ 429566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 439566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match)); 4428b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 45f9720913SMatthew G. Knepley /* Check PLY format */ 469566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 479566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match)); 4828b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 499566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 509566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii)); 519566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig)); 529566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle)); 5328b400f6SJacob Faibussowitsch PetscCheck(!isAscii, PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported"); 54f7d195e4SLawrence Mitchell if (isBinaryLittle) byteSwap = PETSC_TRUE; 559566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 569566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match)); 5728b400f6SJacob 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'; 629566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 639566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match)); 649371c9d4SSatish Balay if (match) 659371c9d4SSatish Balay while (buf != '\n') PetscCall(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR)); 66f9720913SMatthew G. Knepley } 67f9720913SMatthew G. Knepley /* Read vertex information */ 689566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 6928b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 709566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 719566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match)); 7228b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 74f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nv); 7508401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 76f9720913SMatthew G. Knepley match = PETSC_TRUE; 77f9720913SMatthew G. Knepley while (match) { 789566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 799566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match)); 80f9720913SMatthew G. Knepley if (match) { 81f9720913SMatthew G. Knepley PetscBool matchB; 82f9720913SMatthew G. Knepley 839566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 8408401ef6SPierre Jolivet PetscCheck(Nvp < 16, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle more than 16 property statements in PLY file header: %s", line); 85f9720913SMatthew G. Knepley snum = sscanf(line, "%s %s", ntype, name); 8608401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 879566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "float32", 16, &matchB)); 88f9720913SMatthew G. Knepley if (matchB) { 89f9720913SMatthew G. Knepley vtype[Nvp] = 'f'; 90f9720913SMatthew G. Knepley } else { 919566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "int32", 16, &matchB)); 92f9720913SMatthew G. Knepley if (matchB) { 93f9720913SMatthew G. Knepley vtype[Nvp] = 'd'; 94f9720913SMatthew G. Knepley } else { 959566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "uint8", 16, &matchB)); 96f9720913SMatthew G. Knepley if (matchB) { 97f9720913SMatthew G. Knepley vtype[Nvp] = 'c'; 9898921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse type in PLY file header: %s", line); 99f9720913SMatthew G. Knepley } 100f9720913SMatthew G. Knepley } 1019566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "x", 16, &matchB)); 102ad540459SPierre Jolivet if (matchB) xi = Nvp; 1039566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "y", 16, &matchB)); 104ad540459SPierre Jolivet if (matchB) yi = Nvp; 1059566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "z", 16, &matchB)); 106ad540459SPierre Jolivet if (matchB) zi = Nvp; 107f9720913SMatthew G. Knepley ++Nvp; 108f9720913SMatthew G. Knepley } 109f9720913SMatthew G. Knepley } 110f9720913SMatthew G. Knepley /* Read cell information */ 1119566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 11228b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1149566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "face", PETSC_MAX_PATH_LEN, &match)); 11528b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 117f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nc); 11808401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1199566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING)); 120f9720913SMatthew G. Knepley snum = sscanf(line, "property list %s %s %s", ntype, itype, name); 12108401ef6SPierre Jolivet PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1229566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "uint8", 1024, &match)); 12328b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line); 1249566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "vertex_indices", 1024, &match)); 12528b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line); 126f9720913SMatthew G. Knepley /* Header should terminate */ 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1289566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match)); 12928b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file: %s", line); 13025ce1634SJed Brown } else { 13125ce1634SJed Brown Nc = Nv = 0; 132f9720913SMatthew G. Knepley } 1339566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1349566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1359566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, cdim)); 138f9720913SMatthew G. Knepley /* Read coordinates */ 1399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv)); 143f9720913SMatthew G. Knepley for (v = Nc; v < Nc + Nv; ++v) { 1449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, cdim)); 1459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 146f9720913SMatthew G. Knepley } 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1499566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1519566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1529566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, cdim)); 1539566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 1549566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 155dd400576SPatrick Sanan if (rank == 0) { 156f9720913SMatthew G. Knepley float rbuf[1]; 157f9720913SMatthew G. Knepley int ibuf[1]; 158f9720913SMatthew G. Knepley 159f9720913SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 160f9720913SMatthew G. Knepley for (p = 0; p < Nvp; ++p) { 161f9720913SMatthew G. Knepley if (vtype[p] == 'f') { 1629566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT)); 1639566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&rbuf, PETSC_FLOAT, 1)); 164f9720913SMatthew G. Knepley if (p == xi) coords[v * cdim + 0] = rbuf[0]; 165f9720913SMatthew G. Knepley else if (p == yi) coords[v * cdim + 1] = rbuf[0]; 166f9720913SMatthew G. Knepley else if (p == zi) coords[v * cdim + 2] = rbuf[0]; 167f9720913SMatthew G. Knepley } else if (vtype[p] == 'd') { 1689566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT)); 1699566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&ibuf, PETSC_INT, 1)); 170f9720913SMatthew G. Knepley } else if (vtype[p] == 'c') { 1719566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 172f9720913SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid vertex property type in PLY file"); 173f9720913SMatthew G. Knepley } 174f9720913SMatthew G. Knepley } 175f9720913SMatthew G. Knepley } 1769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1779566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 1789566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 179f9720913SMatthew G. Knepley /* Read topology */ 180dd400576SPatrick Sanan if (rank == 0) { 181f9720913SMatthew G. Knepley char ibuf[1]; 182f9720913SMatthew G. Knepley PetscInt vbuf[16], corners; 183f9720913SMatthew G. Knepley 184f9720913SMatthew G. Knepley /* Assume same cells */ 1859566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 186f9720913SMatthew G. Knepley corners = ibuf[0]; 1879566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, corners)); 1889566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 189f9720913SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 19048a46eb9SPierre Jolivet if (c > 0) PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 19163a3b9bcSJacob Faibussowitsch PetscCheck(ibuf[0] == corners, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "All cells must have the same number of vertices in PLY file: %" PetscInt_FMT " != %" PetscInt_FMT, (PetscInt)ibuf[0], corners); 1929566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT)); 1939566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0])); 194f9720913SMatthew G. Knepley for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc; 1959566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, vbuf)); 196f9720913SMatthew G. Knepley } 197f9720913SMatthew G. Knepley } 1989566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 1999566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 201f9720913SMatthew G. Knepley if (interpolate) { 2025fd9971aSMatthew G. Knepley DM idm; 203f9720913SMatthew G. Knepley 2049566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 2059566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 206f9720913SMatthew G. Knepley *dm = idm; 207f9720913SMatthew G. Knepley } 208*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209f9720913SMatthew G. Knepley } 210