xref: /petsc/src/dm/impls/plex/plexply.c (revision 07c2e4feb6773e78bda63e3a89d5b841667f9670)
1f9720913SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
2f9720913SMatthew G. Knepley 
3cc4c1da9SBarry Smith /*@
41d27aa22SBarry Smith   DMPlexCreatePLYFromFile - Create a `DMPLEX` mesh from a PLY <https://en.wikipedia.org/wiki/PLY_(file_format)> file.
5f9720913SMatthew G. Knepley 
6a4e35b19SJacob Faibussowitsch   Input Parameters:
7f9720913SMatthew G. Knepley + comm        - The MPI communicator
86afe31f6SMartin Diehl . filename    - Name of the .ply file
9f9720913SMatthew G. Knepley - interpolate - Create faces and edges in the mesh
10f9720913SMatthew G. Knepley 
11f9720913SMatthew G. Knepley   Output Parameter:
121d27aa22SBarry Smith . dm - The `DMPLEX` object representing the mesh
13f9720913SMatthew G. Knepley 
14f9720913SMatthew G. Knepley   Level: beginner
15f9720913SMatthew G. Knepley 
166afe31f6SMartin Diehl .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()`
17f9720913SMatthew G. Knepley @*/
DMPlexCreatePLYFromFile(MPI_Comm comm,const char filename[],PetscBool interpolate,DM * dm)18d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
19d71ae5a4SJacob Faibussowitsch {
20f9720913SMatthew G. Knepley   PetscViewer  viewer;
21f9720913SMatthew G. Knepley   Vec          coordinates;
22f9720913SMatthew G. Knepley   PetscSection coordSection;
23f9720913SMatthew G. Knepley   PetscScalar *coords;
24f9720913SMatthew G. Knepley   char         line[PETSC_MAX_PATH_LEN], ntype[16], itype[16], name[1024], vtype[16];
25f9720913SMatthew G. Knepley   PetscBool    match, byteSwap = PETSC_FALSE;
2664b21ff4SMatthew G. Knepley   PetscInt     dim = 2, cdim = 3, Nvp = 0, coordSize, xi = -1, yi = -1, zi = -1, v, c, p;
27f9720913SMatthew G. Knepley   PetscMPIInt  rank;
28f9720913SMatthew G. Knepley   int          snum, Nv, Nc;
29f9720913SMatthew G. Knepley 
30f9720913SMatthew G. Knepley   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
329566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, &viewer));
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY));
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(viewer, filename));
36dd400576SPatrick Sanan   if (rank == 0) {
37f9720913SMatthew G. Knepley     PetscBool isAscii, isBinaryBig, isBinaryLittle;
38f9720913SMatthew G. Knepley 
39f9720913SMatthew G. Knepley     /* Check for PLY file */
409566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
419566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match));
4228b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file");
43f9720913SMatthew G. Knepley     /* Check PLY format */
449566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
459566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match));
4628b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file");
479566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
489566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii));
499566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig));
509566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle));
5128b400f6SJacob Faibussowitsch     PetscCheck(!isAscii, PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported");
52f7d195e4SLawrence Mitchell     if (isBinaryLittle) byteSwap = PETSC_TRUE;
539566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
549566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match));
5528b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid version of PLY file, %s", line);
56f9720913SMatthew G. Knepley     /* Ignore comments */
57f9720913SMatthew G. Knepley     match = PETSC_TRUE;
58f9720913SMatthew G. Knepley     while (match) {
59f9720913SMatthew G. Knepley       char buf = '\0';
609566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
619566063dSJacob Faibussowitsch       PetscCall(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match));
629371c9d4SSatish Balay       if (match)
639371c9d4SSatish Balay         while (buf != '\n') PetscCall(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR));
64f9720913SMatthew G. Knepley     }
65f9720913SMatthew G. Knepley     /* Read vertex information */
669566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match));
6728b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
689566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
699566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match));
7028b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
719566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
72f9720913SMatthew G. Knepley     snum = sscanf(line, "%d", &Nv);
7308401ef6SPierre Jolivet     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
74f9720913SMatthew G. Knepley     match = PETSC_TRUE;
75f9720913SMatthew G. Knepley     while (match) {
769566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
779566063dSJacob Faibussowitsch       PetscCall(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match));
78f9720913SMatthew G. Knepley       if (match) {
79f9720913SMatthew G. Knepley         PetscBool matchB;
80f9720913SMatthew G. Knepley 
819566063dSJacob Faibussowitsch         PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING));
8208401ef6SPierre Jolivet         PetscCheck(Nvp < 16, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle more than 16 property statements in PLY file header: %s", line);
83f9720913SMatthew G. Knepley         snum = sscanf(line, "%s %s", ntype, name);
8408401ef6SPierre Jolivet         PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
859566063dSJacob Faibussowitsch         PetscCall(PetscStrncmp(ntype, "float32", 16, &matchB));
86f9720913SMatthew G. Knepley         if (matchB) {
87f9720913SMatthew G. Knepley           vtype[Nvp] = 'f';
88f9720913SMatthew G. Knepley         } else {
899566063dSJacob Faibussowitsch           PetscCall(PetscStrncmp(ntype, "int32", 16, &matchB));
90f9720913SMatthew G. Knepley           if (matchB) {
91f9720913SMatthew G. Knepley             vtype[Nvp] = 'd';
92f9720913SMatthew G. Knepley           } else {
939566063dSJacob Faibussowitsch             PetscCall(PetscStrncmp(ntype, "uint8", 16, &matchB));
94*966bd95aSPierre Jolivet             PetscCheck(matchB, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse type in PLY file header: %s", line);
95f9720913SMatthew G. Knepley             vtype[Nvp] = 'c';
96f9720913SMatthew G. Knepley           }
97f9720913SMatthew G. Knepley         }
989566063dSJacob Faibussowitsch         PetscCall(PetscStrncmp(name, "x", 16, &matchB));
99ad540459SPierre Jolivet         if (matchB) xi = Nvp;
1009566063dSJacob Faibussowitsch         PetscCall(PetscStrncmp(name, "y", 16, &matchB));
101ad540459SPierre Jolivet         if (matchB) yi = Nvp;
1029566063dSJacob Faibussowitsch         PetscCall(PetscStrncmp(name, "z", 16, &matchB));
103ad540459SPierre Jolivet         if (matchB) zi = Nvp;
104f9720913SMatthew G. Knepley         ++Nvp;
105f9720913SMatthew G. Knepley       }
106f9720913SMatthew G. Knepley     }
107f9720913SMatthew G. Knepley     /* Read cell information */
1089566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match));
10928b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
1109566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
1119566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "face", 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));
114f9720913SMatthew G. Knepley     snum = sscanf(line, "%d", &Nc);
11508401ef6SPierre Jolivet     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
1169566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING));
117f9720913SMatthew G. Knepley     snum = sscanf(line, "property list %s %s %s", ntype, itype, name);
11808401ef6SPierre Jolivet     PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
1199566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(ntype, "uint8", 1024, &match));
12028b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line);
1219566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(name, "vertex_indices", 1024, &match));
12228b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line);
123f9720913SMatthew G. Knepley     /* Header should terminate */
1249566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
1259566063dSJacob Faibussowitsch     PetscCall(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match));
12628b400f6SJacob Faibussowitsch     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file: %s", line);
12725ce1634SJed Brown   } else {
12825ce1634SJed Brown     Nc = Nv = 0;
129f9720913SMatthew G. Knepley   }
1309566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
1319566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1329566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv));
1339566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
1349566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDim(*dm, cdim));
135f9720913SMatthew G. Knepley   /* Read coordinates */
1369566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
1379566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
1389566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim));
1399566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv));
140f9720913SMatthew G. Knepley   for (v = Nc; v < Nc + Nv; ++v) {
1419566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, cdim));
1429566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim));
143f9720913SMatthew G. Knepley   }
1449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
1459566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
1469566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
1479566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
1489566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
1499566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, cdim));
1509566063dSJacob Faibussowitsch   PetscCall(VecSetType(coordinates, VECSTANDARD));
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
152dd400576SPatrick Sanan   if (rank == 0) {
153f9720913SMatthew G. Knepley     float rbuf[1];
154f9720913SMatthew G. Knepley     int   ibuf[1];
155f9720913SMatthew G. Knepley 
156f9720913SMatthew G. Knepley     for (v = 0; v < Nv; ++v) {
157f9720913SMatthew G. Knepley       for (p = 0; p < Nvp; ++p) {
158f9720913SMatthew G. Knepley         if (vtype[p] == 'f') {
1599566063dSJacob Faibussowitsch           PetscCall(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT));
1609566063dSJacob Faibussowitsch           if (byteSwap) PetscCall(PetscByteSwap(&rbuf, PETSC_FLOAT, 1));
161f9720913SMatthew G. Knepley           if (p == xi) coords[v * cdim + 0] = rbuf[0];
162f9720913SMatthew G. Knepley           else if (p == yi) coords[v * cdim + 1] = rbuf[0];
163f9720913SMatthew G. Knepley           else if (p == zi) coords[v * cdim + 2] = rbuf[0];
164f9720913SMatthew G. Knepley         } else if (vtype[p] == 'd') {
1659566063dSJacob Faibussowitsch           PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT));
1669566063dSJacob Faibussowitsch           if (byteSwap) PetscCall(PetscByteSwap(&ibuf, PETSC_INT, 1));
167f9720913SMatthew G. Knepley         } else if (vtype[p] == 'c') {
1689566063dSJacob Faibussowitsch           PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
169f9720913SMatthew G. Knepley         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid vertex property type in PLY file");
170f9720913SMatthew G. Knepley       }
171f9720913SMatthew G. Knepley     }
172f9720913SMatthew G. Knepley   }
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
1749566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
1759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
176f9720913SMatthew G. Knepley   /* Read topology */
177dd400576SPatrick Sanan   if (rank == 0) {
178f9720913SMatthew G. Knepley     char     ibuf[1];
179f9720913SMatthew G. Knepley     PetscInt vbuf[16], corners;
180f9720913SMatthew G. Knepley 
181f9720913SMatthew G. Knepley     /* Assume same cells */
1829566063dSJacob Faibussowitsch     PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
183f9720913SMatthew G. Knepley     corners = ibuf[0];
1849566063dSJacob Faibussowitsch     for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, corners));
1859566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dm));
186f9720913SMatthew G. Knepley     for (c = 0; c < Nc; ++c) {
18748a46eb9SPierre Jolivet       if (c > 0) PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
18863a3b9bcSJacob 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);
1899566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT));
1909566063dSJacob Faibussowitsch       if (byteSwap) PetscCall(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0]));
191f9720913SMatthew G. Knepley       for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc;
1929566063dSJacob Faibussowitsch       PetscCall(DMPlexSetCone(*dm, c, vbuf));
193f9720913SMatthew G. Knepley     }
194f9720913SMatthew G. Knepley   }
1959566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
1969566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
1979566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
198f9720913SMatthew G. Knepley   if (interpolate) {
1995fd9971aSMatthew G. Knepley     DM idm;
200f9720913SMatthew G. Knepley 
2019566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
2029566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
203f9720913SMatthew G. Knepley     *dm = idm;
204f9720913SMatthew G. Knepley   }
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206f9720913SMatthew G. Knepley }
207