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