xref: /petsc/src/dm/impls/plex/plexfluent.c (revision 862e4a309d45a165aaa4da0d704ba733429d833a)
1 #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for fileno() */
2 #define PETSCDM_DLL
3 #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
4 
5 /*@C
6   DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file
7 
8 + comm        - The MPI communicator
9 . filename    - Name of the Fluent mesh file
10 - interpolate - Create faces and edges in the mesh
11 
12   Output Parameter:
13 . dm  - The DM object representing the mesh
14 
15   Level: beginner
16 
17 .seealso: `DMPlexCreateFromFile()`, `DMPlexCreateFluent()`, `DMPlexCreate()`
18 @*/
19 PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
20 {
21   PetscViewer viewer;
22 
23   PetscFunctionBegin;
24   /* Create file viewer and build plex */
25   PetscCall(PetscViewerCreate(comm, &viewer));
26   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
27   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
28   PetscCall(PetscViewerFileSetName(viewer, filename));
29   PetscCall(DMPlexCreateFluent(comm, viewer, interpolate, dm));
30   PetscCall(PetscViewerDestroy(&viewer));
31   PetscFunctionReturn(0);
32 }
33 
34 static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim)
35 {
36   PetscInt ret, i = 0;
37 
38   PetscFunctionBegin;
39   do PetscCall(PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR));
40   while (ret > 0 && buffer[i - 1] != '\0' && buffer[i - 1] != delim);
41   if (!ret) buffer[i - 1] = '\0';
42   else buffer[i] = '\0';
43   PetscFunctionReturn(0);
44 }
45 
46 static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary)
47 {
48   int      fdes = 0;
49   FILE    *file;
50   PetscInt i;
51 
52   PetscFunctionBegin;
53   if (binary) {
54     /* Extract raw file descriptor to read binary block */
55     PetscCall(PetscViewerASCIIGetPointer(viewer, &file));
56     fflush(file);
57     fdes = fileno(file);
58   }
59 
60   if (!binary && dtype == PETSC_INT) {
61     char         cbuf[256];
62     unsigned int ibuf;
63     int          snum;
64     /* Parse hexadecimal ascii integers */
65     for (i = 0; i < count; i++) {
66       PetscCall(PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING));
67       snum = sscanf(cbuf, "%x", &ibuf);
68       PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
69       ((PetscInt *)data)[i] = (PetscInt)ibuf;
70     }
71   } else if (binary && dtype == PETSC_INT) {
72     /* Always read 32-bit ints and cast to PetscInt */
73     int *ibuf;
74     PetscCall(PetscMalloc1(count, &ibuf));
75     PetscCall(PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM));
76     PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, count));
77     for (i = 0; i < count; i++) ((PetscInt *)data)[i] = (PetscInt)(ibuf[i]);
78     PetscCall(PetscFree(ibuf));
79 
80   } else if (binary && dtype == PETSC_SCALAR) {
81     float *fbuf;
82     /* Always read 32-bit floats and cast to PetscScalar */
83     PetscCall(PetscMalloc1(count, &fbuf));
84     PetscCall(PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT));
85     PetscCall(PetscByteSwap(fbuf, PETSC_FLOAT, count));
86     for (i = 0; i < count; i++) ((PetscScalar *)data)[i] = (PetscScalar)(fbuf[i]);
87     PetscCall(PetscFree(fbuf));
88   } else {
89     PetscCall(PetscViewerASCIIRead(viewer, data, count, NULL, dtype));
90   }
91   PetscFunctionReturn(0);
92 }
93 
94 static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s)
95 {
96   char buffer[PETSC_MAX_PATH_LEN];
97   int  snum;
98 
99   PetscFunctionBegin;
100   /* Fast-forward to next section and derive its index */
101   PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
102   PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ' '));
103   snum = sscanf(buffer, "%d", &(s->index));
104   /* If we can't match an index return -1 to signal end-of-file */
105   if (snum < 1) {
106     s->index = -1;
107     PetscFunctionReturn(0);
108   }
109 
110   if (s->index == 0) { /* Comment */
111     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
112 
113   } else if (s->index == 2) { /* Dimension */
114     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
115     snum = sscanf(buffer, "%d", &(s->nd));
116     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
117 
118   } else if (s->index == 10 || s->index == 2010) { /* Vertices */
119     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
120     snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
121     PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
122     if (s->zoneID > 0) {
123       PetscInt numCoords = s->last - s->first + 1;
124       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
125       PetscCall(PetscMalloc1(s->nd * numCoords, (PetscScalar **)&s->data));
126       PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd * numCoords, PETSC_SCALAR, s->index == 2010 ? PETSC_TRUE : PETSC_FALSE));
127       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
128     }
129     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
130 
131   } else if (s->index == 12 || s->index == 2012) { /* Cells */
132     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
133     snum = sscanf(buffer, "(%x", &(s->zoneID));
134     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
135     if (s->zoneID == 0) { /* Header section */
136       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
137       PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
138     } else { /* Data section */
139       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
140       PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
141       if (s->nd == 0) {
142         /* Read cell type definitions for mixed cells */
143         PetscInt numCells = s->last - s->first + 1;
144         PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
145         PetscCall(PetscMalloc1(numCells, (PetscInt **)&s->data));
146         PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index == 2012 ? PETSC_TRUE : PETSC_FALSE));
147         PetscCall(PetscFree(s->data));
148         PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
149       }
150     }
151     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
152 
153   } else if (s->index == 13 || s->index == 2013) { /* Faces */
154     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
155     snum = sscanf(buffer, "(%x", &(s->zoneID));
156     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
157     if (s->zoneID == 0) { /* Header section */
158       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
159       PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
160     } else { /* Data section */
161       PetscInt f, numEntries, numFaces;
162       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
163       PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
164       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
165       switch (s->nd) {
166       case 0:
167         numEntries = PETSC_DETERMINE;
168         break;
169       case 2:
170         numEntries = 2 + 2;
171         break; /* linear */
172       case 3:
173         numEntries = 2 + 3;
174         break; /* triangular */
175       case 4:
176         numEntries = 2 + 4;
177         break; /* quadrilateral */
178       default:
179         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file");
180       }
181       numFaces = s->last - s->first + 1;
182       if (numEntries != PETSC_DETERMINE) {
183         /* Allocate space only if we already know the size of the block */
184         PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data));
185       }
186       for (f = 0; f < numFaces; f++) {
187         if (s->nd == 0) {
188           /* Determine the size of the block for "mixed" facets */
189           PetscInt numFaceVert = 0;
190           PetscCall(DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE));
191           if (numEntries == PETSC_DETERMINE) {
192             numEntries = numFaceVert + 2;
193             PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data));
194           } else {
195             PetscCheck(numEntries == numFaceVert + 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files");
196           }
197         }
198         PetscCall(DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt *)s->data)[f * numEntries]), numEntries, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE));
199       }
200       s->nd = numEntries - 2;
201       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
202     }
203     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
204 
205   } else { /* Unknown section type */
206     PetscInt depth = 1;
207     do {
208       /* Match parentheses when parsing unknown sections */
209       do PetscCall(PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR));
210       while (buffer[0] != '(' && buffer[0] != ')');
211       if (buffer[0] == '(') depth++;
212       if (buffer[0] == ')') depth--;
213     } while (depth > 0);
214     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n'));
215   }
216   PetscFunctionReturn(0);
217 }
218 
219 /*@C
220   DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file.
221 
222   Collective
223 
224   Input Parameters:
225 + comm  - The MPI communicator
226 . viewer - The Viewer associated with a Fluent mesh file
227 - interpolate - Create faces and edges in the mesh
228 
229   Output Parameter:
230 . dm  - The DM object representing the mesh
231 
232   Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm
233 
234   Level: beginner
235 
236 .seealso: `DMPLEX`, `DMCreate()`
237 @*/
238 PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
239 {
240   PetscMPIInt  rank;
241   PetscInt     c, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE;
242   PetscInt     numFaces = PETSC_DETERMINE, f, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE;
243   PetscInt    *faces = NULL, *cellVertices = NULL, *faceZoneIDs = NULL;
244   DMLabel      faceSets = NULL;
245   PetscInt     d, coordSize;
246   PetscScalar *coords, *coordsIn = NULL;
247   PetscSection coordSection;
248   Vec          coordinates;
249 
250   PetscFunctionBegin;
251   PetscCallMPI(MPI_Comm_rank(comm, &rank));
252 
253   if (rank == 0) {
254     FluentSection s;
255     do {
256       PetscCall(DMPlexCreateFluent_ReadSection(viewer, &s));
257       if (s.index == 2) { /* Dimension */
258         dim = s.nd;
259 
260       } else if (s.index == 10 || s.index == 2010) { /* Vertices */
261         if (s.zoneID == 0) numVertices = s.last;
262         else {
263           PetscCheck(!coordsIn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files");
264           coordsIn = (PetscScalar *)s.data;
265         }
266 
267       } else if (s.index == 12 || s.index == 2012) { /* Cells */
268         if (s.zoneID == 0) numCells = s.last;
269         else {
270           switch (s.nd) {
271           case 0:
272             numCellVertices = PETSC_DETERMINE;
273             break;
274           case 1:
275             numCellVertices = 3;
276             break; /* triangular */
277           case 2:
278             numCellVertices = 4;
279             break; /* tetrahedral */
280           case 3:
281             numCellVertices = 4;
282             break; /* quadrilateral */
283           case 4:
284             numCellVertices = 8;
285             break; /* hexahedral */
286           case 5:
287             numCellVertices = 5;
288             break; /* pyramid */
289           case 6:
290             numCellVertices = 6;
291             break; /* wedge */
292           default:
293             numCellVertices = PETSC_DETERMINE;
294           }
295         }
296 
297       } else if (s.index == 13 || s.index == 2013) { /* Facets */
298         if (s.zoneID == 0) {                         /* Header section */
299           numFaces = (PetscInt)(s.last - s.first + 1);
300           if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE;
301           else numFaceVertices = s.nd;
302         } else { /* Data section */
303           unsigned int z;
304 
305           PetscCheck(numFaceVertices == PETSC_DETERMINE || s.nd == numFaceVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported");
306           PetscCheck(numFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file");
307           if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd;
308           numFaceEntries = numFaceVertices + 2;
309           if (!faces) PetscCall(PetscMalloc1(numFaces * numFaceEntries, &faces));
310           if (!faceZoneIDs) PetscCall(PetscMalloc1(numFaces, &faceZoneIDs));
311           PetscCall(PetscMemcpy(&faces[(s.first - 1) * numFaceEntries], s.data, (s.last - s.first + 1) * numFaceEntries * sizeof(PetscInt)));
312           /* Record the zoneID for each face set */
313           for (z = s.first - 1; z < s.last; z++) faceZoneIDs[z] = s.zoneID;
314           PetscCall(PetscFree(s.data));
315         }
316       }
317     } while (s.index >= 0);
318   }
319   PetscCallMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm));
320   PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension");
321 
322   /* Allocate cell-vertex mesh */
323   PetscCall(DMCreate(comm, dm));
324   PetscCall(DMSetType(*dm, DMPLEX));
325   PetscCall(DMSetDimension(*dm, dim));
326   PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
327   if (rank == 0) {
328     PetscCheck(numCells >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file");
329     /* If no cell type was given we assume simplices */
330     if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1;
331     for (c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(*dm, c, numCellVertices));
332   }
333   PetscCall(DMSetUp(*dm));
334 
335   if (rank == 0 && faces) {
336     /* Derive cell-vertex list from face-vertex and face-cell maps */
337     PetscCall(PetscMalloc1(numCells * numCellVertices, &cellVertices));
338     for (c = 0; c < numCells * numCellVertices; c++) cellVertices[c] = -1;
339     for (f = 0; f < numFaces; f++) {
340       PetscInt       *cell;
341       const PetscInt  cl   = faces[f * numFaceEntries + numFaceVertices];
342       const PetscInt  cr   = faces[f * numFaceEntries + numFaceVertices + 1];
343       const PetscInt *face = &(faces[f * numFaceEntries]);
344 
345       if (cl > 0) {
346         cell = &(cellVertices[(cl - 1) * numCellVertices]);
347         for (v = 0; v < numFaceVertices; v++) {
348           PetscBool found = PETSC_FALSE;
349           for (c = 0; c < numCellVertices; c++) {
350             if (cell[c] < 0) break;
351             if (cell[c] == face[v] - 1 + numCells) {
352               found = PETSC_TRUE;
353               break;
354             }
355           }
356           if (!found) cell[c] = face[v] - 1 + numCells;
357         }
358       }
359       if (cr > 0) {
360         cell = &(cellVertices[(cr - 1) * numCellVertices]);
361         for (v = 0; v < numFaceVertices; v++) {
362           PetscBool found = PETSC_FALSE;
363           for (c = 0; c < numCellVertices; c++) {
364             if (cell[c] < 0) break;
365             if (cell[c] == face[v] - 1 + numCells) {
366               found = PETSC_TRUE;
367               break;
368             }
369           }
370           if (!found) cell[c] = face[v] - 1 + numCells;
371         }
372       }
373     }
374     for (c = 0; c < numCells; c++) PetscCall(DMPlexSetCone(*dm, c, &(cellVertices[c * numCellVertices])));
375   }
376   PetscCall(DMPlexSymmetrize(*dm));
377   PetscCall(DMPlexStratify(*dm));
378   if (interpolate) {
379     DM idm;
380 
381     PetscCall(DMPlexInterpolate(*dm, &idm));
382     PetscCall(DMDestroy(dm));
383     *dm = idm;
384   }
385 
386   if (rank == 0 && faces) {
387     PetscInt        fi, joinSize, meetSize, *fverts, cells[2];
388     const PetscInt *join, *meet;
389     PetscCall(PetscMalloc1(numFaceVertices, &fverts));
390     /* Mark facets by finding the full join of all adjacent vertices */
391     for (f = 0; f < numFaces; f++) {
392       const PetscInt cl = faces[f * numFaceEntries + numFaceVertices] - 1;
393       const PetscInt cr = faces[f * numFaceEntries + numFaceVertices + 1] - 1;
394       if (cl > 0 && cr > 0) {
395         /* If we know both adjoining cells we can use a single-level meet */
396         cells[0] = cl;
397         cells[1] = cr;
398         PetscCall(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet));
399         PetscCheck(meetSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f);
400         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], faceZoneIDs[f]));
401         PetscCall(DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet));
402       } else {
403         for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f * numFaceEntries + fi] + numCells - 1;
404         PetscCall(DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join));
405         PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f);
406         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], faceZoneIDs[f]));
407         PetscCall(DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join));
408       }
409     }
410     PetscCall(PetscFree(fverts));
411   }
412 
413   { /* Create Face Sets label at all processes */
414     enum {
415       n = 1
416     };
417     PetscBool flag[n];
418 
419     flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE;
420     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
421     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
422   }
423 
424   /* Read coordinates */
425   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
426   PetscCall(PetscSectionSetNumFields(coordSection, 1));
427   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
428   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
429   for (v = numCells; v < numCells + numVertices; ++v) {
430     PetscCall(PetscSectionSetDof(coordSection, v, dim));
431     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
432   }
433   PetscCall(PetscSectionSetUp(coordSection));
434   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
435   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
436   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
437   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
438   PetscCall(VecSetType(coordinates, VECSTANDARD));
439   PetscCall(VecGetArray(coordinates, &coords));
440   if (rank == 0 && coordsIn) {
441     for (v = 0; v < numVertices; ++v) {
442       for (d = 0; d < dim; ++d) coords[v * dim + d] = coordsIn[v * dim + d];
443     }
444   }
445   PetscCall(VecRestoreArray(coordinates, &coords));
446   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
447   PetscCall(VecDestroy(&coordinates));
448 
449   if (rank == 0) {
450     PetscCall(PetscFree(cellVertices));
451     PetscCall(PetscFree(faces));
452     PetscCall(PetscFree(faceZoneIDs));
453     PetscCall(PetscFree(coordsIn));
454   }
455   PetscFunctionReturn(0);
456 }
457