Actual source code: plexply.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>

  4: /*@C
  5:   DMPlexCreatePLYFromFile - Create a `DMPLEX` mesh from a PLY <https://en.wikipedia.org/wiki/PLY_(file_format)> file.

  7:   Input Parameters:
  8: + comm        - The MPI communicator
  9: . filename    - Name of the .med file
 10: - interpolate - Create faces and edges in the mesh

 12:   Output Parameter:
 13: . dm - The `DMPLEX` object representing the mesh

 15:   Level: beginner

 17: .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateMedFromFile()`, `DMPlexCreateGmsh()`, `DMPlexCreate()`
 18: @*/
 19: PetscErrorCode DMPlexCreatePLYFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
 20: {
 21:   PetscViewer  viewer;
 22:   Vec          coordinates;
 23:   PetscSection coordSection;
 24:   PetscScalar *coords;
 25:   char         line[PETSC_MAX_PATH_LEN], ntype[16], itype[16], name[1024], vtype[16];
 26:   PetscBool    match, byteSwap = PETSC_FALSE;
 27:   PetscInt     dim = 2, cdim = 3, Nvp = 0, coordSize, xi = -1, yi = -1, zi = -1, v, c, p;
 28:   PetscMPIInt  rank;
 29:   int          snum, Nv, Nc;

 31:   PetscFunctionBegin;
 32:   PetscCallMPI(MPI_Comm_rank(comm, &rank));
 33:   PetscCall(PetscViewerCreate(comm, &viewer));
 34:   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY));
 35:   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
 36:   PetscCall(PetscViewerFileSetName(viewer, filename));
 37:   if (rank == 0) {
 38:     PetscBool isAscii, isBinaryBig, isBinaryLittle;

 40:     /* Check for PLY file */
 41:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 42:     PetscCall(PetscStrncmp(line, "ply", PETSC_MAX_PATH_LEN, &match));
 43:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file");
 44:     /* Check PLY format */
 45:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 46:     PetscCall(PetscStrncmp(line, "format", PETSC_MAX_PATH_LEN, &match));
 47:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file");
 48:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 49:     PetscCall(PetscStrncmp(line, "ascii", PETSC_MAX_PATH_LEN, &isAscii));
 50:     PetscCall(PetscStrncmp(line, "binary_big_endian", PETSC_MAX_PATH_LEN, &isBinaryBig));
 51:     PetscCall(PetscStrncmp(line, "binary_little_endian", PETSC_MAX_PATH_LEN, &isBinaryLittle));
 52:     PetscCheck(!isAscii, PETSC_COMM_SELF, PETSC_ERR_SUP, "PLY ascii format not yet supported");
 53:     if (isBinaryLittle) byteSwap = PETSC_TRUE;
 54:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 55:     PetscCall(PetscStrncmp(line, "1.0", PETSC_MAX_PATH_LEN, &match));
 56:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid version of PLY file, %s", line);
 57:     /* Ignore comments */
 58:     match = PETSC_TRUE;
 59:     while (match) {
 60:       char buf = '\0';
 61:       PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 62:       PetscCall(PetscStrncmp(line, "comment", PETSC_MAX_PATH_LEN, &match));
 63:       if (match)
 64:         while (buf != '\n') PetscCall(PetscViewerRead(viewer, &buf, 1, NULL, PETSC_CHAR));
 65:     }
 66:     /* Read vertex information */
 67:     PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match));
 68:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
 69:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 70:     PetscCall(PetscStrncmp(line, "vertex", PETSC_MAX_PATH_LEN, &match));
 71:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
 72:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 73:     snum = sscanf(line, "%d", &Nv);
 74:     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
 75:     match = PETSC_TRUE;
 76:     while (match) {
 77:       PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
 78:       PetscCall(PetscStrncmp(line, "property", PETSC_MAX_PATH_LEN, &match));
 79:       if (match) {
 80:         PetscBool matchB;

 82:         PetscCall(PetscViewerRead(viewer, line, 2, NULL, PETSC_STRING));
 83:         PetscCheck(Nvp < 16, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cannot handle more than 16 property statements in PLY file header: %s", line);
 84:         snum = sscanf(line, "%s %s", ntype, name);
 85:         PetscCheck(snum == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
 86:         PetscCall(PetscStrncmp(ntype, "float32", 16, &matchB));
 87:         if (matchB) {
 88:           vtype[Nvp] = 'f';
 89:         } else {
 90:           PetscCall(PetscStrncmp(ntype, "int32", 16, &matchB));
 91:           if (matchB) {
 92:             vtype[Nvp] = 'd';
 93:           } else {
 94:             PetscCall(PetscStrncmp(ntype, "uint8", 16, &matchB));
 95:             if (matchB) {
 96:               vtype[Nvp] = 'c';
 97:             } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse type in PLY file header: %s", line);
 98:           }
 99:         }
100:         PetscCall(PetscStrncmp(name, "x", 16, &matchB));
101:         if (matchB) xi = Nvp;
102:         PetscCall(PetscStrncmp(name, "y", 16, &matchB));
103:         if (matchB) yi = Nvp;
104:         PetscCall(PetscStrncmp(name, "z", 16, &matchB));
105:         if (matchB) zi = Nvp;
106:         ++Nvp;
107:       }
108:     }
109:     /* Read cell information */
110:     PetscCall(PetscStrncmp(line, "element", PETSC_MAX_PATH_LEN, &match));
111:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
112:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
113:     PetscCall(PetscStrncmp(line, "face", PETSC_MAX_PATH_LEN, &match));
114:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
115:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
116:     snum = sscanf(line, "%d", &Nc);
117:     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
118:     PetscCall(PetscViewerRead(viewer, line, 5, NULL, PETSC_STRING));
119:     snum = sscanf(line, "property list %s %s %s", ntype, itype, name);
120:     PetscCheck(snum == 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to parse PLY file header: %s", line);
121:     PetscCall(PetscStrncmp(ntype, "uint8", 1024, &match));
122:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid size type in PLY file header: %s", line);
123:     PetscCall(PetscStrncmp(name, "vertex_indices", 1024, &match));
124:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid property in PLY file header: %s", line);
125:     /* Header should terminate */
126:     PetscCall(PetscViewerRead(viewer, line, 1, NULL, PETSC_STRING));
127:     PetscCall(PetscStrncmp(line, "end_header", PETSC_MAX_PATH_LEN, &match));
128:     PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid PLY file: %s", line);
129:   } else {
130:     Nc = Nv = 0;
131:   }
132:   PetscCall(DMCreate(comm, dm));
133:   PetscCall(DMSetType(*dm, DMPLEX));
134:   PetscCall(DMPlexSetChart(*dm, 0, Nc + Nv));
135:   PetscCall(DMSetDimension(*dm, dim));
136:   PetscCall(DMSetCoordinateDim(*dm, cdim));
137:   /* Read coordinates */
138:   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
139:   PetscCall(PetscSectionSetNumFields(coordSection, 1));
140:   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, cdim));
141:   PetscCall(PetscSectionSetChart(coordSection, Nc, Nc + Nv));
142:   for (v = Nc; v < Nc + Nv; ++v) {
143:     PetscCall(PetscSectionSetDof(coordSection, v, cdim));
144:     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, cdim));
145:   }
146:   PetscCall(PetscSectionSetUp(coordSection));
147:   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
148:   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
149:   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
150:   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
151:   PetscCall(VecSetBlockSize(coordinates, cdim));
152:   PetscCall(VecSetType(coordinates, VECSTANDARD));
153:   PetscCall(VecGetArray(coordinates, &coords));
154:   if (rank == 0) {
155:     float rbuf[1];
156:     int   ibuf[1];

158:     for (v = 0; v < Nv; ++v) {
159:       for (p = 0; p < Nvp; ++p) {
160:         if (vtype[p] == 'f') {
161:           PetscCall(PetscViewerRead(viewer, &rbuf, 1, NULL, PETSC_FLOAT));
162:           if (byteSwap) PetscCall(PetscByteSwap(&rbuf, PETSC_FLOAT, 1));
163:           if (p == xi) coords[v * cdim + 0] = rbuf[0];
164:           else if (p == yi) coords[v * cdim + 1] = rbuf[0];
165:           else if (p == zi) coords[v * cdim + 2] = rbuf[0];
166:         } else if (vtype[p] == 'd') {
167:           PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_INT));
168:           if (byteSwap) PetscCall(PetscByteSwap(&ibuf, PETSC_INT, 1));
169:         } else if (vtype[p] == 'c') {
170:           PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
171:         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid vertex property type in PLY file");
172:       }
173:     }
174:   }
175:   PetscCall(VecRestoreArray(coordinates, &coords));
176:   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
177:   PetscCall(VecDestroy(&coordinates));
178:   /* Read topology */
179:   if (rank == 0) {
180:     char     ibuf[1];
181:     PetscInt vbuf[16], corners;

183:     /* Assume same cells */
184:     PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
185:     corners = ibuf[0];
186:     for (c = 0; c < Nc; ++c) PetscCall(DMPlexSetConeSize(*dm, c, corners));
187:     PetscCall(DMSetUp(*dm));
188:     for (c = 0; c < Nc; ++c) {
189:       if (c > 0) PetscCall(PetscViewerRead(viewer, &ibuf, 1, NULL, PETSC_CHAR));
190:       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);
191:       PetscCall(PetscViewerRead(viewer, &vbuf, ibuf[0], NULL, PETSC_INT));
192:       if (byteSwap) PetscCall(PetscByteSwap(&vbuf, PETSC_INT, ibuf[0]));
193:       for (v = 0; v < ibuf[0]; ++v) vbuf[v] += Nc;
194:       PetscCall(DMPlexSetCone(*dm, c, vbuf));
195:     }
196:   }
197:   PetscCall(DMPlexSymmetrize(*dm));
198:   PetscCall(DMPlexStratify(*dm));
199:   PetscCall(PetscViewerDestroy(&viewer));
200:   if (interpolate) {
201:     DM idm;

203:     PetscCall(DMPlexInterpolate(*dm, &idm));
204:     PetscCall(DMDestroy(dm));
205:     *dm = idm;
206:   }
207:   PetscFunctionReturn(PETSC_SUCCESS);
208: }