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 5*1d27aa22SBarry Smith DMPlexCreatePLYFromFile - Create a `DMPLEX` mesh from a PLY <https://en.wikipedia.org/wiki/PLY_(file_format)> file. 6f9720913SMatthew G. Knepley 7a4e35b19SJacob Faibussowitsch Input Parameters: 8f9720913SMatthew G. Knepley + comm - The MPI communicator 9f9720913SMatthew G. Knepley . filename - Name of the .med file 10f9720913SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 11f9720913SMatthew G. Knepley 12f9720913SMatthew G. Knepley Output Parameter: 13*1d27aa22SBarry Smith . dm - The `DMPLEX` object representing the mesh 14f9720913SMatthew G. Knepley 15f9720913SMatthew G. Knepley Level: beginner 16f9720913SMatthew G. Knepley 17db781477SPatrick Sanan .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 18f9720913SMatthew G. Knepley @*/ 19d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 20d71ae5a4SJacob Faibussowitsch { 21f9720913SMatthew G. Knepley PetscViewer viewer; 22f9720913SMatthew G. Knepley Vec coordinates; 23f9720913SMatthew G. Knepley PetscSection coordSection; 24f9720913SMatthew G. Knepley PetscScalar *coords; 25f9720913SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN], ntype[16], itype[16], name[1024], vtype[16]; 26f9720913SMatthew G. Knepley PetscBool match, byteSwap = PETSC_FALSE; 2764b21ff4SMatthew G. Knepley PetscInt dim = 2, cdim = 3, Nvp = 0, coordSize, xi = -1, yi = -1, zi = -1, v, c, p; 28f9720913SMatthew G. Knepley PetscMPIInt rank; 29f9720913SMatthew G. Knepley int snum, Nv, Nc; 30f9720913SMatthew G. Knepley 31f9720913SMatthew G. Knepley PetscFunctionBegin; 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY)); 359566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 37dd400576SPatrick Sanan if (rank == 0) { 38f9720913SMatthew G. Knepley PetscBool isAscii, isBinaryBig, isBinaryLittle; 39f9720913SMatthew G. Knepley 40f9720913SMatthew G. Knepley /* Check for PLY file */ 419566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 429566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match)); 4328b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 44f9720913SMatthew G. Knepley /* Check PLY format */ 459566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 469566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match)); 4728b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 499566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii)); 509566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig)); 519566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle)); 5228b400f6SJacob Faibussowitsch PetscCheck(!isAscii, PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported"); 53f7d195e4SLawrence Mitchell if (isBinaryLittle) byteSwap = PETSC_TRUE; 549566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 559566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match)); 5628b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid version of PLY file, %s", line); 57f9720913SMatthew G. Knepley /* Ignore comments */ 58f9720913SMatthew G. Knepley match = PETSC_TRUE; 59f9720913SMatthew G. Knepley while (match) { 60f9720913SMatthew G. Knepley char buf = '\0'; 619566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 629566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match)); 639371c9d4SSatish Balay if (match) 649371c9d4SSatish Balay while (buf != '\n') PetscCall(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR)); 65f9720913SMatthew G. Knepley } 66f9720913SMatthew G. Knepley /* Read vertex information */ 679566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 6828b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 709566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match)); 7128b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 73f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nv); 7408401ef6SPierre Jolivet PetscCheck(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) { 779566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 789566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match)); 79f9720913SMatthew G. Knepley if (match) { 80f9720913SMatthew G. Knepley PetscBool matchB; 81f9720913SMatthew G. Knepley 829566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 8308401ef6SPierre Jolivet PetscCheck(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); 8508401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 869566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "float32", 16, &matchB)); 87f9720913SMatthew G. Knepley if (matchB) { 88f9720913SMatthew G. Knepley vtype[Nvp] = 'f'; 89f9720913SMatthew G. Knepley } else { 909566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "int32", 16, &matchB)); 91f9720913SMatthew G. Knepley if (matchB) { 92f9720913SMatthew G. Knepley vtype[Nvp] = 'd'; 93f9720913SMatthew G. Knepley } else { 949566063dSJacob Faibussowitsch PetscCall(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 } 1009566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "x", 16, &matchB)); 101ad540459SPierre Jolivet if (matchB) xi = Nvp; 1029566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "y", 16, &matchB)); 103ad540459SPierre Jolivet if (matchB) yi = Nvp; 1049566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "z", 16, &matchB)); 105ad540459SPierre Jolivet if (matchB) zi = Nvp; 106f9720913SMatthew G. Knepley ++Nvp; 107f9720913SMatthew G. Knepley } 108f9720913SMatthew G. Knepley } 109f9720913SMatthew G. Knepley /* Read cell information */ 1109566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match)); 11128b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1139566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "face", PETSC_MAX_PATH_LEN, &match)); 11428b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 116f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nc); 11708401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1189566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING)); 119f9720913SMatthew G. Knepley snum = sscanf(line, "property list %s %s %s", ntype, itype, name); 12008401ef6SPierre Jolivet PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1219566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "uint8", 1024, &match)); 12228b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line); 1239566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "vertex_indices", 1024, &match)); 12428b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line); 125f9720913SMatthew G. Knepley /* Header should terminate */ 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1279566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match)); 12828b400f6SJacob 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 } 1329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1349566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv)); 1359566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, cdim)); 137f9720913SMatthew G. Knepley /* Read coordinates */ 1389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 1399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 1419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv)); 142f9720913SMatthew G. Knepley for (v = Nc; v < Nc + Nv; ++v) { 1439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, cdim)); 1449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 145f9720913SMatthew G. Knepley } 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1489566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1519566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, cdim)); 1529566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 1539566063dSJacob Faibussowitsch PetscCall(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') { 1619566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT)); 1629566063dSJacob Faibussowitsch if (byteSwap) PetscCall(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') { 1679566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT)); 1689566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&ibuf, PETSC_INT, 1)); 169f9720913SMatthew G. Knepley } else if (vtype[p] == 'c') { 1709566063dSJacob Faibussowitsch PetscCall(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 } 1759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1769566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 1779566063dSJacob Faibussowitsch PetscCall(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 */ 1849566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 185f9720913SMatthew G. Knepley corners = ibuf[0]; 1869566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, corners)); 1879566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 188f9720913SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 18948a46eb9SPierre Jolivet if (c > 0) PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 19063a3b9bcSJacob 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); 1919566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT)); 1929566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0])); 193f9720913SMatthew G. Knepley for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc; 1949566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, vbuf)); 195f9720913SMatthew G. Knepley } 196f9720913SMatthew G. Knepley } 1979566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 1989566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 200f9720913SMatthew G. Knepley if (interpolate) { 2015fd9971aSMatthew G. Knepley DM idm; 202f9720913SMatthew G. Knepley 2039566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 2049566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 205f9720913SMatthew G. Knepley *dm = idm; 206f9720913SMatthew G. Knepley } 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 208f9720913SMatthew G. Knepley } 209