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*a4e35b19SJacob Faibussowitsch DMPlexCreatePLYFromFile - Create a `DMPLEX` mesh from a PLY file. 6f9720913SMatthew G. Knepley 7*a4e35b19SJacob 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: 13f9720913SMatthew G. Knepley . dm - The DM object representing the mesh 14f9720913SMatthew G. Knepley 15*a4e35b19SJacob Faibussowitsch References: 16*a4e35b19SJacob Faibussowitsch . * - https://en.wikipedia.org/wiki/PLY_(file_format) 17f9720913SMatthew G. Knepley 18f9720913SMatthew G. Knepley Level: beginner 19f9720913SMatthew G. Knepley 20db781477SPatrick Sanan .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()` 21f9720913SMatthew G. Knepley @*/ 22d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 23d71ae5a4SJacob Faibussowitsch { 24f9720913SMatthew G. Knepley PetscViewer viewer; 25f9720913SMatthew G. Knepley Vec coordinates; 26f9720913SMatthew G. Knepley PetscSection coordSection; 27f9720913SMatthew G. Knepley PetscScalar *coords; 28f9720913SMatthew G. Knepley char line[PETSC_MAX_PATH_LEN], ntype[16], itype[16], name[1024], vtype[16]; 29f9720913SMatthew G. Knepley PetscBool match, byteSwap = PETSC_FALSE; 3064b21ff4SMatthew G. Knepley PetscInt dim = 2, cdim = 3, Nvp = 0, coordSize, xi = -1, yi = -1, zi = -1, v, c, p; 31f9720913SMatthew G. Knepley PetscMPIInt rank; 32f9720913SMatthew G. Knepley int snum, Nv, Nc; 33f9720913SMatthew G. Knepley 34f9720913SMatthew G. Knepley PetscFunctionBegin; 359566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 369566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 379566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 399566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 40dd400576SPatrick Sanan if (rank == 0) { 41f9720913SMatthew G. Knepley PetscBool isAscii, isBinaryBig, isBinaryLittle; 42f9720913SMatthew G. Knepley 43f9720913SMatthew G. Knepley /* Check for PLY file */ 449566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 459566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match)); 4628b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 47f9720913SMatthew G. Knepley /* Check PLY format */ 489566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 499566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match)); 5028b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file"); 519566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 529566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii)); 539566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig)); 549566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle)); 5528b400f6SJacob Faibussowitsch PetscCheck(!isAscii, PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported"); 56f7d195e4SLawrence Mitchell if (isBinaryLittle) byteSwap = PETSC_TRUE; 579566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 589566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match)); 5928b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid version of PLY file, %s", line); 60f9720913SMatthew G. Knepley /* Ignore comments */ 61f9720913SMatthew G. Knepley match = PETSC_TRUE; 62f9720913SMatthew G. Knepley while (match) { 63f9720913SMatthew G. Knepley char buf = '\0'; 649566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 659566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match)); 669371c9d4SSatish Balay if (match) 679371c9d4SSatish Balay while (buf != '\n') PetscCall(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR)); 68f9720913SMatthew G. Knepley } 69f9720913SMatthew G. Knepley /* Read vertex information */ 709566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", 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)); 739566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match)); 7428b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 76f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nv); 7708401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 78f9720913SMatthew G. Knepley match = PETSC_TRUE; 79f9720913SMatthew G. Knepley while (match) { 809566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 819566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match)); 82f9720913SMatthew G. Knepley if (match) { 83f9720913SMatthew G. Knepley PetscBool matchB; 84f9720913SMatthew G. Knepley 859566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING)); 8608401ef6SPierre Jolivet PetscCheck(Nvp < 16, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle more than 16 property statements in PLY file header: %s", line); 87f9720913SMatthew G. Knepley snum = sscanf(line, "%s %s", ntype, name); 8808401ef6SPierre Jolivet PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 899566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "float32", 16, &matchB)); 90f9720913SMatthew G. Knepley if (matchB) { 91f9720913SMatthew G. Knepley vtype[Nvp] = 'f'; 92f9720913SMatthew G. Knepley } else { 939566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "int32", 16, &matchB)); 94f9720913SMatthew G. Knepley if (matchB) { 95f9720913SMatthew G. Knepley vtype[Nvp] = 'd'; 96f9720913SMatthew G. Knepley } else { 979566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "uint8", 16, &matchB)); 98f9720913SMatthew G. Knepley if (matchB) { 99f9720913SMatthew G. Knepley vtype[Nvp] = 'c'; 10098921bdaSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse type in PLY file header: %s", line); 101f9720913SMatthew G. Knepley } 102f9720913SMatthew G. Knepley } 1039566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "x", 16, &matchB)); 104ad540459SPierre Jolivet if (matchB) xi = Nvp; 1059566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "y", 16, &matchB)); 106ad540459SPierre Jolivet if (matchB) yi = Nvp; 1079566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "z", 16, &matchB)); 108ad540459SPierre Jolivet if (matchB) zi = Nvp; 109f9720913SMatthew G. Knepley ++Nvp; 110f9720913SMatthew G. Knepley } 111f9720913SMatthew G. Knepley } 112f9720913SMatthew G. Knepley /* Read cell information */ 1139566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "element", 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)); 1169566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "face", PETSC_MAX_PATH_LEN, &match)); 11728b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1189566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 119f9720913SMatthew G. Knepley snum = sscanf(line, "%d", &Nc); 12008401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING)); 122f9720913SMatthew G. Knepley snum = sscanf(line, "property list %s %s %s", ntype, itype, name); 12308401ef6SPierre Jolivet PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line); 1249566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(ntype, "uint8", 1024, &match)); 12528b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line); 1269566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(name, "vertex_indices", 1024, &match)); 12728b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line); 128f9720913SMatthew G. Knepley /* Header should terminate */ 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING)); 1309566063dSJacob Faibussowitsch PetscCall(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match)); 13128b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file: %s", line); 13225ce1634SJed Brown } else { 13325ce1634SJed Brown Nc = Nv = 0; 134f9720913SMatthew G. Knepley } 1359566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 1369566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1379566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv)); 1389566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 1399566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, cdim)); 140f9720913SMatthew G. Knepley /* Read coordinates */ 1419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim)); 1449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv)); 145f9720913SMatthew G. Knepley for (v = Nc; v < Nc + Nv; ++v) { 1469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, cdim)); 1479566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim)); 148f9720913SMatthew G. Knepley } 1499566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 1509566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1519566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1529566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1539566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1549566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, cdim)); 1559566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 1569566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 157dd400576SPatrick Sanan if (rank == 0) { 158f9720913SMatthew G. Knepley float rbuf[1]; 159f9720913SMatthew G. Knepley int ibuf[1]; 160f9720913SMatthew G. Knepley 161f9720913SMatthew G. Knepley for (v = 0; v < Nv; ++v) { 162f9720913SMatthew G. Knepley for (p = 0; p < Nvp; ++p) { 163f9720913SMatthew G. Knepley if (vtype[p] == 'f') { 1649566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT)); 1659566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&rbuf, PETSC_FLOAT, 1)); 166f9720913SMatthew G. Knepley if (p == xi) coords[v * cdim + 0] = rbuf[0]; 167f9720913SMatthew G. Knepley else if (p == yi) coords[v * cdim + 1] = rbuf[0]; 168f9720913SMatthew G. Knepley else if (p == zi) coords[v * cdim + 2] = rbuf[0]; 169f9720913SMatthew G. Knepley } else if (vtype[p] == 'd') { 1709566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT)); 1719566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&ibuf, PETSC_INT, 1)); 172f9720913SMatthew G. Knepley } else if (vtype[p] == 'c') { 1739566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 174f9720913SMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid vertex property type in PLY file"); 175f9720913SMatthew G. Knepley } 176f9720913SMatthew G. Knepley } 177f9720913SMatthew G. Knepley } 1789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 1799566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 181f9720913SMatthew G. Knepley /* Read topology */ 182dd400576SPatrick Sanan if (rank == 0) { 183f9720913SMatthew G. Knepley char ibuf[1]; 184f9720913SMatthew G. Knepley PetscInt vbuf[16], corners; 185f9720913SMatthew G. Knepley 186f9720913SMatthew G. Knepley /* Assume same cells */ 1879566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 188f9720913SMatthew G. Knepley corners = ibuf[0]; 1899566063dSJacob Faibussowitsch for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, corners)); 1909566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 191f9720913SMatthew G. Knepley for (c = 0; c < Nc; ++c) { 19248a46eb9SPierre Jolivet if (c > 0) PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR)); 19363a3b9bcSJacob 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); 1949566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT)); 1959566063dSJacob Faibussowitsch if (byteSwap) PetscCall(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0])); 196f9720913SMatthew G. Knepley for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc; 1979566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, vbuf)); 198f9720913SMatthew G. Knepley } 199f9720913SMatthew G. Knepley } 2009566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 2019566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 203f9720913SMatthew G. Knepley if (interpolate) { 2045fd9971aSMatthew G. Knepley DM idm; 205f9720913SMatthew G. Knepley 2069566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 2079566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 208f9720913SMatthew G. Knepley *dm = idm; 209f9720913SMatthew G. Knepley } 2103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 211f9720913SMatthew G. Knepley } 212