xref: /petsc/src/dm/impls/plex/plexfluent.c (revision 07c2e4feb6773e78bda63e3a89d5b841667f9670)
10039db0dSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for fileno() */
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
32f0bd6dcSMichael Lange 
45733454dSMatthew G. Knepley /* Utility struct to store the contents of a Fluent file in memory */
55733454dSMatthew G. Knepley typedef struct {
65733454dSMatthew G. Knepley   int          index; /* Type of section */
75733454dSMatthew G. Knepley   unsigned int zoneID;
85733454dSMatthew G. Knepley   unsigned int first;
95733454dSMatthew G. Knepley   unsigned int last;
105733454dSMatthew G. Knepley   int          type;
115733454dSMatthew G. Knepley   int          nd; /* Either ND or element-type */
125733454dSMatthew G. Knepley   void        *data;
135733454dSMatthew G. Knepley } FluentSection;
145733454dSMatthew G. Knepley 
15cc4c1da9SBarry Smith /*@
16a1cb98faSBarry Smith   DMPlexCreateFluentFromFile - Create a `DMPLEX` mesh from a Fluent mesh file
172f0bd6dcSMichael Lange 
18a1cb98faSBarry Smith   Collective
19a1cb98faSBarry Smith 
20a1cb98faSBarry Smith   Input Parameters:
212f0bd6dcSMichael Lange + comm        - The MPI communicator
222f0bd6dcSMichael Lange . filename    - Name of the Fluent mesh file
232f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh
242f0bd6dcSMichael Lange 
252f0bd6dcSMichael Lange   Output Parameter:
26a1cb98faSBarry Smith . dm - The `DM` object representing the mesh
272f0bd6dcSMichael Lange 
282f0bd6dcSMichael Lange   Level: beginner
292f0bd6dcSMichael Lange 
301cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateFluent()`, `DMPlexCreate()`
312f0bd6dcSMichael Lange @*/
DMPlexCreateFluentFromFile(MPI_Comm comm,const char filename[],PetscBool interpolate,DM * dm)32d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
33d71ae5a4SJacob Faibussowitsch {
342f0bd6dcSMichael Lange   PetscViewer viewer;
352f0bd6dcSMichael Lange 
362f0bd6dcSMichael Lange   PetscFunctionBegin;
372f0bd6dcSMichael Lange   /* Create file viewer and build plex */
389566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, &viewer));
399566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
409566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
419566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(viewer, filename));
429566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFluent(comm, viewer, interpolate, dm));
439566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
452f0bd6dcSMichael Lange }
462f0bd6dcSMichael Lange 
DMPlexCreateFluent_ReadString(PetscViewer viewer,char * buffer,char delim)47d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim)
48d71ae5a4SJacob Faibussowitsch {
494b375a39SMichael Lange   PetscInt ret, i = 0;
502f0bd6dcSMichael Lange 
512f0bd6dcSMichael Lange   PetscFunctionBegin;
52f4f49eeaSPierre Jolivet   do PetscCall(PetscViewerRead(viewer, &buffer[i++], 1, &ret, PETSC_CHAR));
53c7cc6ee6SMatthew G. Knepley   while (ret > 0 && buffer[i - 1] != '\0' && buffer[i - 1] != delim && i < PETSC_MAX_PATH_LEN - 1);
549371c9d4SSatish Balay   if (!ret) buffer[i - 1] = '\0';
559371c9d4SSatish Balay   else buffer[i] = '\0';
56c7cc6ee6SMatthew G. Knepley   PetscCheck(i < PETSC_MAX_PATH_LEN - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Buffer overflow! This is not a valid Fluent file.");
573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
582f0bd6dcSMichael Lange }
592f0bd6dcSMichael Lange 
DMPlexCreateFluent_ReadValues(PetscViewer viewer,void * data,PetscInt count,PetscDataType dtype,PetscBool binary,PetscInt * numClosingParens)60c7cc6ee6SMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary, PetscInt *numClosingParens)
61d71ae5a4SJacob Faibussowitsch {
6270c9a859SSatish Balay   int      fdes = 0;
6319d58f9dSMichael Lange   FILE    *file;
6419d58f9dSMichael Lange   PetscInt i;
6519d58f9dSMichael Lange 
6619d58f9dSMichael Lange   PetscFunctionBegin;
67c7cc6ee6SMatthew G. Knepley   *numClosingParens = 0;
6819d58f9dSMichael Lange   if (binary) {
6919d58f9dSMichael Lange     /* Extract raw file descriptor to read binary block */
709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetPointer(viewer, &file));
71c69effb2SJacob Faibussowitsch     PetscCall(PetscFFlush(file));
729371c9d4SSatish Balay     fdes = fileno(file);
7319d58f9dSMichael Lange   }
7419d58f9dSMichael Lange 
7519d58f9dSMichael Lange   if (!binary && dtype == PETSC_INT) {
7619d58f9dSMichael Lange     char         cbuf[256];
77cfb60857SMatthew G. Knepley     unsigned int ibuf;
78cfb60857SMatthew G. Knepley     int          snum;
7919d58f9dSMichael Lange     /* Parse hexadecimal ascii integers */
8019d58f9dSMichael Lange     for (i = 0; i < count; i++) {
81c7cc6ee6SMatthew G. Knepley       size_t len;
82c7cc6ee6SMatthew G. Knepley 
839566063dSJacob Faibussowitsch       PetscCall(PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING));
8419d58f9dSMichael Lange       snum = sscanf(cbuf, "%x", &ibuf);
8508401ef6SPierre Jolivet       PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
8619d58f9dSMichael Lange       ((PetscInt *)data)[i] = (PetscInt)ibuf;
87c7cc6ee6SMatthew G. Knepley       // Check for trailing parentheses
88c7cc6ee6SMatthew G. Knepley       PetscCall(PetscStrlen(cbuf, &len));
89c7cc6ee6SMatthew G. Knepley       while (cbuf[len - 1] == ')' && len > 0) {
90c7cc6ee6SMatthew G. Knepley         ++(*numClosingParens);
91c7cc6ee6SMatthew G. Knepley         --len;
92c7cc6ee6SMatthew G. Knepley       }
9319d58f9dSMichael Lange     }
9419d58f9dSMichael Lange   } else if (binary && dtype == PETSC_INT) {
9519d58f9dSMichael Lange     /* Always read 32-bit ints and cast to PetscInt */
9619d58f9dSMichael Lange     int *ibuf;
979566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, &ibuf));
989566063dSJacob Faibussowitsch     PetscCall(PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM));
999566063dSJacob Faibussowitsch     PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, count));
100835f2295SStefano Zampini     for (i = 0; i < count; i++) ((PetscInt *)data)[i] = ibuf[i];
1019566063dSJacob Faibussowitsch     PetscCall(PetscFree(ibuf));
10219d58f9dSMichael Lange 
10319d58f9dSMichael Lange   } else if (binary && dtype == PETSC_SCALAR) {
10419d58f9dSMichael Lange     float *fbuf;
10519d58f9dSMichael Lange     /* Always read 32-bit floats and cast to PetscScalar */
1069566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(count, &fbuf));
1079566063dSJacob Faibussowitsch     PetscCall(PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT));
1089566063dSJacob Faibussowitsch     PetscCall(PetscByteSwap(fbuf, PETSC_FLOAT, count));
109835f2295SStefano Zampini     for (i = 0; i < count; i++) ((PetscScalar *)data)[i] = fbuf[i];
1109566063dSJacob Faibussowitsch     PetscCall(PetscFree(fbuf));
11119d58f9dSMichael Lange   } else {
1129566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIRead(viewer, data, count, NULL, dtype));
11319d58f9dSMichael Lange   }
1143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11519d58f9dSMichael Lange }
11619d58f9dSMichael Lange 
DMPlexCreateFluent_ReadSection(PetscViewer viewer,FluentSection * s)117d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s)
118d71ae5a4SJacob Faibussowitsch {
1192f0bd6dcSMichael Lange   char buffer[PETSC_MAX_PATH_LEN];
1202f0bd6dcSMichael Lange   int  snum;
1212f0bd6dcSMichael Lange 
1222f0bd6dcSMichael Lange   PetscFunctionBegin;
1232f0bd6dcSMichael Lange   /* Fast-forward to next section and derive its index */
1249566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
1259566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ' '));
126f4f49eeaSPierre Jolivet   snum = sscanf(buffer, "%d", &s->index);
1272f0bd6dcSMichael Lange   /* If we can't match an index return -1 to signal end-of-file */
1289371c9d4SSatish Balay   if (snum < 1) {
1299371c9d4SSatish Balay     s->index = -1;
1303ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
1319371c9d4SSatish Balay   }
1322f0bd6dcSMichael Lange 
1332f0bd6dcSMichael Lange   if (s->index == 0) { /* Comment */
1349566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
1352f0bd6dcSMichael Lange 
1362f0bd6dcSMichael Lange   } else if (s->index == 2) { /* Dimension */
1379566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
138f4f49eeaSPierre Jolivet     snum = sscanf(buffer, "%d", &s->nd);
13908401ef6SPierre Jolivet     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1402f0bd6dcSMichael Lange 
14119d58f9dSMichael Lange   } else if (s->index == 10 || s->index == 2010) { /* Vertices */
142c7cc6ee6SMatthew G. Knepley     PetscInt numClosingParens = 0;
143c7cc6ee6SMatthew G. Knepley 
1449566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
145f4f49eeaSPierre Jolivet     snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd);
14608401ef6SPierre Jolivet     PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1473f6dc66eSMichael Lange     if (s->zoneID > 0) {
1483f6dc66eSMichael Lange       PetscInt numCoords = s->last - s->first + 1;
1499566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
1509566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(s->nd * numCoords, (PetscScalar **)&s->data));
151c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd * numCoords, PETSC_SCALAR, s->index == 2010 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens));
152c7cc6ee6SMatthew G. Knepley       if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
153c7cc6ee6SMatthew G. Knepley       else --numClosingParens;
1543f6dc66eSMichael Lange     }
155c7cc6ee6SMatthew G. Knepley     if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
156c7cc6ee6SMatthew G. Knepley     else --numClosingParens;
157c7cc6ee6SMatthew G. Knepley     PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
15819d58f9dSMichael Lange   } else if (s->index == 12 || s->index == 2012) { /* Cells */
159c7cc6ee6SMatthew G. Knepley     PetscInt numClosingParens = 0;
160c7cc6ee6SMatthew G. Knepley 
1619566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
162f4f49eeaSPierre Jolivet     snum = sscanf(buffer, "(%x", &s->zoneID);
16308401ef6SPierre Jolivet     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1642f0bd6dcSMichael Lange     if (s->zoneID == 0) { /* Header section */
165f4f49eeaSPierre Jolivet       snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd);
16608401ef6SPierre Jolivet       PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
167f7320561SMichael Lange     } else { /* Data section */
168f4f49eeaSPierre Jolivet       snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd);
16908401ef6SPierre Jolivet       PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
170895fcc3eSMichael Lange       if (s->nd == 0) {
171895fcc3eSMichael Lange         /* Read cell type definitions for mixed cells */
17219d58f9dSMichael Lange         PetscInt numCells = s->last - s->first + 1;
1739566063dSJacob Faibussowitsch         PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
1749566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(numCells, (PetscInt **)&s->data));
175c7cc6ee6SMatthew G. Knepley         PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index == 2012 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens));
176c7cc6ee6SMatthew G. Knepley         if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
177c7cc6ee6SMatthew G. Knepley         else --numClosingParens;
178895fcc3eSMichael Lange       }
1792f0bd6dcSMichael Lange     }
180c7cc6ee6SMatthew G. Knepley     if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
181c7cc6ee6SMatthew G. Knepley     else --numClosingParens;
182c7cc6ee6SMatthew G. Knepley     PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
18319d58f9dSMichael Lange   } else if (s->index == 13 || s->index == 2013) { /* Faces */
184c7cc6ee6SMatthew G. Knepley     PetscInt numClosingParens = 0;
185c7cc6ee6SMatthew G. Knepley 
1869566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
187f4f49eeaSPierre Jolivet     snum = sscanf(buffer, "(%x", &s->zoneID);
18808401ef6SPierre Jolivet     PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1892f0bd6dcSMichael Lange     if (s->zoneID == 0) { /* Header section */
190f4f49eeaSPierre Jolivet       snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd);
19108401ef6SPierre Jolivet       PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
192f7320561SMichael Lange     } else { /* Data section */
1935733454dSMatthew G. Knepley       PetscInt numEntries, numFaces, maxsize = 0, offset = 0;
1945733454dSMatthew G. Knepley 
195f4f49eeaSPierre Jolivet       snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd);
19608401ef6SPierre Jolivet       PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1979566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '('));
198f7320561SMichael Lange       switch (s->nd) {
199d71ae5a4SJacob Faibussowitsch       case 0:
200d71ae5a4SJacob Faibussowitsch         numEntries = PETSC_DETERMINE;
201d71ae5a4SJacob Faibussowitsch         break;
202d71ae5a4SJacob Faibussowitsch       case 2:
203d71ae5a4SJacob Faibussowitsch         numEntries = 2 + 2;
204d71ae5a4SJacob Faibussowitsch         break; /* linear */
205d71ae5a4SJacob Faibussowitsch       case 3:
206d71ae5a4SJacob Faibussowitsch         numEntries = 2 + 3;
207d71ae5a4SJacob Faibussowitsch         break; /* triangular */
208d71ae5a4SJacob Faibussowitsch       case 4:
209d71ae5a4SJacob Faibussowitsch         numEntries = 2 + 4;
210d71ae5a4SJacob Faibussowitsch         break; /* quadrilateral */
211d71ae5a4SJacob Faibussowitsch       default:
212d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file");
213f7320561SMichael Lange       }
214f7320561SMichael Lange       numFaces = s->last - s->first + 1;
215895fcc3eSMichael Lange       if (numEntries != PETSC_DETERMINE) {
216895fcc3eSMichael Lange         /* Allocate space only if we already know the size of the block */
2179566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data));
218895fcc3eSMichael Lange       }
2195733454dSMatthew G. Knepley       for (PetscInt f = 0; f < numFaces; f++) {
220895fcc3eSMichael Lange         if (s->nd == 0) {
221895fcc3eSMichael Lange           /* Determine the size of the block for "mixed" facets */
222137d0469SJed Brown           PetscInt numFaceVert = 0;
223c7cc6ee6SMatthew G. Knepley           PetscCall(DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens));
2245733454dSMatthew G. Knepley           if (!f) {
2255733454dSMatthew G. Knepley             maxsize = (numFaceVert + 3) * numFaces;
2265733454dSMatthew G. Knepley             PetscCall(PetscMalloc1(maxsize, (PetscInt **)&s->data));
227895fcc3eSMichael Lange           } else {
2285733454dSMatthew G. Knepley             if (offset + numFaceVert + 3 >= maxsize) {
2295733454dSMatthew G. Knepley               PetscInt *tmp;
2305733454dSMatthew G. Knepley 
2315733454dSMatthew G. Knepley               PetscCall(PetscMalloc1(maxsize * 2, &tmp));
2325733454dSMatthew G. Knepley               PetscCall(PetscArraycpy(tmp, (PetscInt *)s->data, maxsize));
2335733454dSMatthew G. Knepley               PetscCall(PetscFree(s->data));
2345733454dSMatthew G. Knepley               maxsize *= 2;
2355733454dSMatthew G. Knepley               s->data = tmp;
236895fcc3eSMichael Lange             }
237895fcc3eSMichael Lange           }
2385733454dSMatthew G. Knepley           ((PetscInt *)s->data)[offset] = numFaceVert;
2395733454dSMatthew G. Knepley           ++offset;
2405733454dSMatthew G. Knepley           numEntries = numFaceVert + 2;
241f7320561SMichael Lange         }
2425733454dSMatthew G. Knepley         PetscCall(DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt *)s->data)[offset]), numEntries, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens));
2435733454dSMatthew G. Knepley         offset += numEntries;
2445733454dSMatthew G. Knepley       }
2455733454dSMatthew G. Knepley       if (s->nd != 0) PetscCall(PetscMPIIntCast(numEntries - 2, &s->nd));
246c7cc6ee6SMatthew G. Knepley       if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
247c7cc6ee6SMatthew G. Knepley       else --numClosingParens;
2482f0bd6dcSMichael Lange     }
249c7cc6ee6SMatthew G. Knepley     if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
250c7cc6ee6SMatthew G. Knepley     else --numClosingParens;
251c7cc6ee6SMatthew G. Knepley     PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
252c7cc6ee6SMatthew G. Knepley   } else if (s->index == 39) { /* Label information */
253c7cc6ee6SMatthew G. Knepley     char labelName[PETSC_MAX_PATH_LEN];
254c7cc6ee6SMatthew G. Knepley     char caseName[PETSC_MAX_PATH_LEN];
2552f0bd6dcSMichael Lange 
256c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')'));
257c3279546SPierre Jolivet     snum = sscanf(buffer, "(%u %s %s %d)", &s->zoneID, caseName, labelName, &s->nd);
258c7cc6ee6SMatthew G. Knepley     PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file: %d", snum);
259c7cc6ee6SMatthew G. Knepley     PetscInt depth = 1;
260c7cc6ee6SMatthew G. Knepley     do {
261c7cc6ee6SMatthew G. Knepley       /* Match parentheses when parsing unknown sections */
262c7cc6ee6SMatthew G. Knepley       do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR));
263c7cc6ee6SMatthew G. Knepley       while (buffer[0] != '(' && buffer[0] != ')');
264c7cc6ee6SMatthew G. Knepley       if (buffer[0] == '(') depth++;
265c7cc6ee6SMatthew G. Knepley       if (buffer[0] == ')') depth--;
266c7cc6ee6SMatthew G. Knepley     } while (depth > 0);
267c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n'));
268c7cc6ee6SMatthew G. Knepley     PetscCall(PetscStrallocpy(labelName, (char **)&s->data));
269c3279546SPierre Jolivet     PetscCall(PetscInfo((PetscObject)viewer, "CASE: Zone ID %u is label %s\n", s->zoneID, labelName));
2702f0bd6dcSMichael Lange   } else { /* Unknown section type */
271060da220SMatthew G. Knepley     PetscInt depth = 1;
2722f0bd6dcSMichael Lange     do {
2732f0bd6dcSMichael Lange       /* Match parentheses when parsing unknown sections */
274f4f49eeaSPierre Jolivet       do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR));
2752f0bd6dcSMichael Lange       while (buffer[0] != '(' && buffer[0] != ')');
2762f0bd6dcSMichael Lange       if (buffer[0] == '(') depth++;
2772f0bd6dcSMichael Lange       if (buffer[0] == ')') depth--;
2782f0bd6dcSMichael Lange     } while (depth > 0);
2799566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n'));
2802f0bd6dcSMichael Lange   }
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2822f0bd6dcSMichael Lange }
2832f0bd6dcSMichael Lange 
284c7cc6ee6SMatthew G. Knepley // Inserts point `face` with orientation `ornt` into the cone of point `cell` at position `c`, which is the first empty slot
InsertFace(DM dm,PetscInt cell,PetscInt face,PetscInt ornt)285c7cc6ee6SMatthew G. Knepley static PetscErrorCode InsertFace(DM dm, PetscInt cell, PetscInt face, PetscInt ornt)
286c7cc6ee6SMatthew G. Knepley {
287c7cc6ee6SMatthew G. Knepley   const PetscInt *cone;
288c7cc6ee6SMatthew G. Knepley   PetscInt        coneSize, c;
289c7cc6ee6SMatthew G. Knepley 
290c7cc6ee6SMatthew G. Knepley   PetscFunctionBegin;
291c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetCone(dm, cell, &cone));
292c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
293c7cc6ee6SMatthew G. Knepley   for (c = 0; c < coneSize; ++c)
294c7cc6ee6SMatthew G. Knepley     if (cone[c] < 0) break;
2955733454dSMatthew G. Knepley   PetscCheck(c < coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be inserted in cone of cell %" PetscInt_FMT " with size %" PetscInt_FMT, face, cell, coneSize);
296c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexInsertCone(dm, cell, c, face));
297c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexInsertConeOrientation(dm, cell, c, ornt));
298c7cc6ee6SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
299c7cc6ee6SMatthew G. Knepley }
300c7cc6ee6SMatthew G. Knepley 
ReorderPolygon(DM dm,PetscInt cell)301c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderPolygon(DM dm, PetscInt cell)
302c7cc6ee6SMatthew G. Knepley {
303c7cc6ee6SMatthew G. Knepley   const PetscInt *cone, *ornt;
304c7cc6ee6SMatthew G. Knepley   PetscInt        coneSize, newCone[16], newOrnt[16];
305c7cc6ee6SMatthew G. Knepley 
306c7cc6ee6SMatthew G. Knepley   PetscFunctionBegin;
307c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt));
308c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetConeSize(dm, cell, &coneSize));
309c7cc6ee6SMatthew G. Knepley   newCone[0] = cone[0];
310c7cc6ee6SMatthew G. Knepley   newOrnt[0] = ornt[0];
311c7cc6ee6SMatthew G. Knepley   for (PetscInt c = 1; c < coneSize; ++c) {
312c7cc6ee6SMatthew G. Knepley     const PetscInt *fcone;
313c7cc6ee6SMatthew G. Knepley     PetscInt        firstVertex, lastVertex, c2;
314c7cc6ee6SMatthew G. Knepley 
315c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexGetCone(dm, newCone[c - 1], &fcone));
316c7cc6ee6SMatthew G. Knepley     lastVertex = newOrnt[c - 1] ? fcone[0] : fcone[1];
317c7cc6ee6SMatthew G. Knepley     for (c2 = 0; c2 < coneSize; ++c2) {
318c7cc6ee6SMatthew G. Knepley       const PetscInt *fcone2;
319c7cc6ee6SMatthew G. Knepley 
320c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexGetCone(dm, cone[c2], &fcone2));
321c7cc6ee6SMatthew G. Knepley       firstVertex = ornt[c2] ? fcone2[1] : fcone2[0];
322c7cc6ee6SMatthew G. Knepley       if (lastVertex == firstVertex) {
323c7cc6ee6SMatthew G. Knepley         // Point `cell` matched point `lastVertex` on face `cone[c2]` with orientation `ornt[c2]`
324c7cc6ee6SMatthew G. Knepley         break;
325c7cc6ee6SMatthew G. Knepley       }
326c7cc6ee6SMatthew G. Knepley     }
327c7cc6ee6SMatthew G. Knepley     PetscCheck(c2 < coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match as position %" PetscInt_FMT, cell, c);
328c7cc6ee6SMatthew G. Knepley     newCone[c] = cone[c2];
329c7cc6ee6SMatthew G. Knepley     newOrnt[c] = ornt[c2];
330c7cc6ee6SMatthew G. Knepley   }
331c7cc6ee6SMatthew G. Knepley   {
332c7cc6ee6SMatthew G. Knepley     const PetscInt *fcone, *fcone2;
333c7cc6ee6SMatthew G. Knepley     PetscInt        vertex, vertex2;
334c7cc6ee6SMatthew G. Knepley 
335c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexGetCone(dm, newCone[coneSize - 1], &fcone));
336c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexGetCone(dm, newCone[0], &fcone2));
337c7cc6ee6SMatthew G. Knepley     vertex  = newOrnt[coneSize - 1] ? fcone[0] : fcone[1];
338c7cc6ee6SMatthew G. Knepley     vertex2 = newOrnt[0] ? fcone2[1] : fcone2[0];
339c7cc6ee6SMatthew G. Knepley     PetscCheck(vertex == vertex2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " did not match at the endpoint", cell);
340c7cc6ee6SMatthew G. Knepley   }
341c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetCone(dm, cell, newCone));
342c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt));
343c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt));
344c7cc6ee6SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
345c7cc6ee6SMatthew G. Knepley }
346c7cc6ee6SMatthew G. Knepley 
ReorderTetrahedron(PetscViewer viewer,DM dm,PetscInt cell)347fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderTetrahedron(PetscViewer viewer, DM dm, PetscInt cell)
348c7cc6ee6SMatthew G. Knepley {
349c7cc6ee6SMatthew G. Knepley   const PetscInt *cone, *ornt, *fcone, *fornt, *farr, faces[4] = {0, 1, 3, 2};
350c7cc6ee6SMatthew G. Knepley   PetscInt        newCone[16], newOrnt[16];
351c7cc6ee6SMatthew G. Knepley 
352c7cc6ee6SMatthew G. Knepley   PetscFunctionBegin;
353c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt));
354c7cc6ee6SMatthew G. Knepley   newCone[0] = cone[0];
355c7cc6ee6SMatthew G. Knepley   newOrnt[0] = ornt[0];
356c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt));
357c7cc6ee6SMatthew G. Knepley   farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]);
358c7cc6ee6SMatthew G. Knepley   // Loop over each edge in the initial triangle
359c7cc6ee6SMatthew G. Knepley   for (PetscInt e = 0; e < 3; ++e) {
360c7cc6ee6SMatthew G. Knepley     const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]);
361c7cc6ee6SMatthew G. Knepley     PetscInt       c;
362c7cc6ee6SMatthew G. Knepley 
363c7cc6ee6SMatthew G. Knepley     // Loop over each remaining face in the tetrahedron
364c7cc6ee6SMatthew G. Knepley     //   On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face
365c7cc6ee6SMatthew G. Knepley     for (c = 1; c < 4; ++c) {
366c7cc6ee6SMatthew G. Knepley       const PetscInt *fcone2, *fornt2, *farr2;
367c7cc6ee6SMatthew G. Knepley       PetscInt        c2;
368fc8f8d65SMatthew G. Knepley       PetscBool       flip = PETSC_FALSE;
369c7cc6ee6SMatthew G. Knepley 
370c7cc6ee6SMatthew G. Knepley       // Checking face `cone[c]` with orientation `ornt[c]`
371c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2));
372c7cc6ee6SMatthew G. Knepley       farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]);
373c7cc6ee6SMatthew G. Knepley       // Check for edge
374c7cc6ee6SMatthew G. Knepley       for (c2 = 0; c2 < 3; ++c2) {
375fc8f8d65SMatthew G. Knepley         const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]);
376c7cc6ee6SMatthew G. Knepley         // Trying to match edge `edge2` with final orientation `eornt2`
377c7cc6ee6SMatthew G. Knepley         if (edge == edge2) {
378fc8f8d65SMatthew G. Knepley           //PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell % " PetscInt_FMT " edge %" PetscInt_FMT " (%" PetscInt_FMT ") found twice with the same orientation in face %" PetscInt_FMT " edge %" PetscInt_FMT, cell, edge, e, c, c2);
379c7cc6ee6SMatthew G. Knepley           // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge`
380fc8f8d65SMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Matched cell %" PetscInt_FMT " edge %" PetscInt_FMT "/%" PetscInt_FMT " (%" PetscInt_FMT ") to face %" PetscInt_FMT "/%" PetscInt_FMT " edge %" PetscInt_FMT " (%" PetscInt_FMT ")\n", cell, edge, e, eornt, cone[c], c, c2, eornt2));
381fc8f8d65SMatthew G. Knepley           flip = eornt != -(eornt2 + 1) ? PETSC_TRUE : PETSC_FALSE;
382c7cc6ee6SMatthew G. Knepley           break;
383c7cc6ee6SMatthew G. Knepley         }
384c7cc6ee6SMatthew G. Knepley       }
385c7cc6ee6SMatthew G. Knepley       if (c2 < 3) {
386c7cc6ee6SMatthew G. Knepley         newCone[faces[e + 1]] = cone[c];
387c7cc6ee6SMatthew G. Knepley         // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0)
388fc8f8d65SMatthew G. Knepley         //   Face 1 should match its edge 2
389fc8f8d65SMatthew G. Knepley         //   Face 2 should match its edge 0
390fc8f8d65SMatthew G. Knepley         //   Face 3 should match its edge 0
391fc8f8d65SMatthew G. Knepley         if (flip) {
392fc8f8d65SMatthew G. Knepley           newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, -((c2 + (!e ? 1 : 2)) % 3 + 1), ornt[c]);
393fc8f8d65SMatthew G. Knepley         } else {
394c7cc6ee6SMatthew G. Knepley           newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, !e ? (c2 + 1) % 3 : c2, ornt[c]);
395fc8f8d65SMatthew G. Knepley         }
396c7cc6ee6SMatthew G. Knepley         break;
397c7cc6ee6SMatthew G. Knepley       }
398c7cc6ee6SMatthew G. Knepley     }
399c7cc6ee6SMatthew G. Knepley     PetscCheck(c < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e);
400c7cc6ee6SMatthew G. Knepley   }
401c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt));
402c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetCone(dm, cell, newCone));
403c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt));
404c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt));
405c7cc6ee6SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
406c7cc6ee6SMatthew G. Knepley }
407c7cc6ee6SMatthew G. Knepley 
ReorderHexahedron(DM dm,PetscInt cell)408c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderHexahedron(DM dm, PetscInt cell)
409c7cc6ee6SMatthew G. Knepley {
410c7cc6ee6SMatthew G. Knepley   const PetscInt *cone, *ornt, *fcone, *fornt, *farr;
411c7cc6ee6SMatthew G. Knepley   const PetscInt  faces[6] = {0, 5, 3, 4, 2, 1};
412c7cc6ee6SMatthew G. Knepley   PetscInt        used[6]  = {1, 0, 0, 0, 0, 0};
413c7cc6ee6SMatthew G. Knepley   PetscInt        newCone[16], newOrnt[16];
414c7cc6ee6SMatthew G. Knepley 
415c7cc6ee6SMatthew G. Knepley   PetscFunctionBegin;
416c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt));
417c7cc6ee6SMatthew G. Knepley   newCone[0] = cone[0];
418c7cc6ee6SMatthew G. Knepley   newOrnt[0] = ornt[0];
419c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt));
420c7cc6ee6SMatthew G. Knepley   farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[0]);
421c7cc6ee6SMatthew G. Knepley   // Loop over each edge in the initial quadrilateral
422c7cc6ee6SMatthew G. Knepley   for (PetscInt e = 0; e < 4; ++e) {
423c7cc6ee6SMatthew G. Knepley     const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]);
424c7cc6ee6SMatthew G. Knepley     PetscInt       c;
425c7cc6ee6SMatthew G. Knepley 
426c7cc6ee6SMatthew G. Knepley     // Loop over each remaining face in the hexahedron
427c7cc6ee6SMatthew G. Knepley     //   On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face
428c7cc6ee6SMatthew G. Knepley     for (c = 1; c < 6; ++c) {
429c7cc6ee6SMatthew G. Knepley       const PetscInt *fcone2, *fornt2, *farr2;
430c7cc6ee6SMatthew G. Knepley       PetscInt        c2;
431c7cc6ee6SMatthew G. Knepley 
432c7cc6ee6SMatthew G. Knepley       // Checking face `cone[c]` with orientation `ornt[c]`
433c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2));
434c7cc6ee6SMatthew G. Knepley       farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]);
435c7cc6ee6SMatthew G. Knepley       // Check for edge
436c7cc6ee6SMatthew G. Knepley       for (c2 = 0; c2 < 4; ++c2) {
437c7cc6ee6SMatthew G. Knepley         const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]);
438c7cc6ee6SMatthew G. Knepley         // Trying to match edge `edge2` with final orientation `eornt2`
439c7cc6ee6SMatthew G. Knepley         if (edge == edge2) {
440c7cc6ee6SMatthew G. Knepley           PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge);
441c7cc6ee6SMatthew G. Knepley           // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge`
442c7cc6ee6SMatthew G. Knepley           break;
443c7cc6ee6SMatthew G. Knepley         }
444c7cc6ee6SMatthew G. Knepley       }
445c7cc6ee6SMatthew G. Knepley       if (c2 < 4) {
446c7cc6ee6SMatthew G. Knepley         used[c]               = 1;
447c7cc6ee6SMatthew G. Knepley         newCone[faces[e + 1]] = cone[c];
448c7cc6ee6SMatthew G. Knepley         // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0)
449c7cc6ee6SMatthew G. Knepley         newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, !e ? (c2 + 1) % 4 : c2, ornt[c]);
450c7cc6ee6SMatthew G. Knepley         break;
451c7cc6ee6SMatthew G. Knepley       }
452c7cc6ee6SMatthew G. Knepley     }
453c7cc6ee6SMatthew G. Knepley     PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e);
454c7cc6ee6SMatthew G. Knepley   }
455c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt));
456c7cc6ee6SMatthew G. Knepley   // Add last face
457c7cc6ee6SMatthew G. Knepley   {
458c7cc6ee6SMatthew G. Knepley     PetscInt c, c2;
459c7cc6ee6SMatthew G. Knepley 
460c7cc6ee6SMatthew G. Knepley     for (c = 1; c < 6; ++c)
461c7cc6ee6SMatthew G. Knepley       if (!used[c]) break;
462c7cc6ee6SMatthew G. Knepley     PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell);
463c7cc6ee6SMatthew G. Knepley     // Match first edge to 3rd edge in newCone[2]
464c7cc6ee6SMatthew G. Knepley     {
465c7cc6ee6SMatthew G. Knepley       const PetscInt *fcone2, *fornt2, *farr2;
466c7cc6ee6SMatthew G. Knepley 
467c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt));
468c7cc6ee6SMatthew G. Knepley       farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]);
469c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2));
470c7cc6ee6SMatthew G. Knepley       farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]);
471c7cc6ee6SMatthew G. Knepley 
472c7cc6ee6SMatthew G. Knepley       const PetscInt e    = 2;
473c7cc6ee6SMatthew G. Knepley       const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]);
474c7cc6ee6SMatthew G. Knepley       // Trying to match edge `edge` with final orientation `eornt` of face `newCone[2]` to some edge of face `cone[c]` with orientation `ornt[c]`
475c7cc6ee6SMatthew G. Knepley       for (c2 = 0; c2 < 4; ++c2) {
476c7cc6ee6SMatthew G. Knepley         const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]);
477c7cc6ee6SMatthew G. Knepley         // Trying to match edge `edge2` with final orientation `eornt2`
478c7cc6ee6SMatthew G. Knepley         if (edge == edge2) {
479c7cc6ee6SMatthew G. Knepley           PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge);
480c7cc6ee6SMatthew G. Knepley           // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge`
481c7cc6ee6SMatthew G. Knepley           break;
482c7cc6ee6SMatthew G. Knepley         }
483c7cc6ee6SMatthew G. Knepley       }
484c7cc6ee6SMatthew G. Knepley       PetscCheck(c2 < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in");
485c7cc6ee6SMatthew G. Knepley     }
486c7cc6ee6SMatthew G. Knepley     newCone[faces[5]] = cone[c];
487c7cc6ee6SMatthew G. Knepley     // Compute new orientation of face based on which edge was matched
488c7cc6ee6SMatthew G. Knepley     newOrnt[faces[5]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]);
489c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt));
490c7cc6ee6SMatthew G. Knepley   }
491c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetCone(dm, cell, newCone));
492c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt));
493c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt));
494c7cc6ee6SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
495c7cc6ee6SMatthew G. Knepley }
496c7cc6ee6SMatthew G. Knepley 
4975733454dSMatthew G. Knepley // {0, 1, 2}, {3, 4, 5}, {0, 2, 4, 3}, {2, 1, 5, 4}, {1, 0, 3, 5}
ReorderWedge(DM dm,PetscInt cell)4985733454dSMatthew G. Knepley static PetscErrorCode ReorderWedge(DM dm, PetscInt cell)
4995733454dSMatthew G. Knepley {
5005733454dSMatthew G. Knepley   const PetscInt *cone, *ornt, *fcone, *fornt, *farr;
5015733454dSMatthew G. Knepley   const PetscInt  faces[5] = {0, 4, 3, 2, 1};
5025733454dSMatthew G. Knepley   PetscInt        used[5]  = {0, 0, 0, 0, 0};
5035733454dSMatthew G. Knepley   PetscInt        newCone[16], newOrnt[16], cS, bottom = 0;
5045733454dSMatthew G. Knepley 
5055733454dSMatthew G. Knepley   PetscFunctionBegin;
5065733454dSMatthew G. Knepley   PetscCall(DMPlexGetConeSize(dm, cell, &cS));
5075733454dSMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt));
5085733454dSMatthew G. Knepley   for (PetscInt c = 0; c < cS; ++c) {
5095733454dSMatthew G. Knepley     DMPolytopeType ct;
5105733454dSMatthew G. Knepley 
5115733454dSMatthew G. Knepley     PetscCall(DMPlexGetCellType(dm, cone[c], &ct));
5125733454dSMatthew G. Knepley     if (ct == DM_POLYTOPE_TRIANGLE) {
5135733454dSMatthew G. Knepley       bottom = c;
5145733454dSMatthew G. Knepley       break;
5155733454dSMatthew G. Knepley     }
5165733454dSMatthew G. Knepley   }
5175733454dSMatthew G. Knepley   used[bottom] = 1;
5185733454dSMatthew G. Knepley   newCone[0]   = cone[bottom];
5195733454dSMatthew G. Knepley   newOrnt[0]   = ornt[bottom];
5205733454dSMatthew G. Knepley   PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt));
5215733454dSMatthew G. Knepley   farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]);
5225733454dSMatthew G. Knepley   // Loop over each edge in the initial triangle
5235733454dSMatthew G. Knepley   for (PetscInt e = 0; e < 3; ++e) {
5245733454dSMatthew G. Knepley     const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]);
5255733454dSMatthew G. Knepley     PetscInt       c;
5265733454dSMatthew G. Knepley 
5275733454dSMatthew G. Knepley     // Loop over each remaining face in the prism
5285733454dSMatthew G. Knepley     //   On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face
5295733454dSMatthew G. Knepley     for (c = 0; c < 5; ++c) {
5305733454dSMatthew G. Knepley       const PetscInt *fcone2, *fornt2, *farr2;
5315733454dSMatthew G. Knepley       DMPolytopeType  ct;
5325733454dSMatthew G. Knepley       PetscInt        c2;
5335733454dSMatthew G. Knepley 
5345733454dSMatthew G. Knepley       if (c == bottom) continue;
5355733454dSMatthew G. Knepley       PetscCall(DMPlexGetCellType(dm, cone[c], &ct));
5365733454dSMatthew G. Knepley       if (ct != DM_POLYTOPE_QUADRILATERAL) continue;
5375733454dSMatthew G. Knepley       // Checking face `cone[c]` with orientation `ornt[c]`
5385733454dSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2));
5395733454dSMatthew G. Knepley       farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]);
5405733454dSMatthew G. Knepley       // Check for edge
5415733454dSMatthew G. Knepley       for (c2 = 0; c2 < 4; ++c2) {
5425733454dSMatthew G. Knepley         const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]);
5435733454dSMatthew G. Knepley         // Trying to match edge `edge2` with final orientation `eornt2`
5445733454dSMatthew G. Knepley         if (edge == edge2) {
5455733454dSMatthew G. Knepley           PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge);
5465733454dSMatthew G. Knepley           // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge`
5475733454dSMatthew G. Knepley           break;
5485733454dSMatthew G. Knepley         }
5495733454dSMatthew G. Knepley       }
5505733454dSMatthew G. Knepley       if (c2 < 4) {
5515733454dSMatthew G. Knepley         used[c]               = 1;
5525733454dSMatthew G. Knepley         newCone[faces[e + 1]] = cone[c];
5535733454dSMatthew G. Knepley         // Compute new orientation of face based on which edge was matched, edge 0 should always match the bottom
5545733454dSMatthew G. Knepley         newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]);
5555733454dSMatthew G. Knepley         break;
5565733454dSMatthew G. Knepley       }
5575733454dSMatthew G. Knepley     }
5585733454dSMatthew G. Knepley     PetscCheck(c < 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e);
5595733454dSMatthew G. Knepley   }
5605733454dSMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt));
5615733454dSMatthew G. Knepley   // Add last face
5625733454dSMatthew G. Knepley   {
5635733454dSMatthew G. Knepley     PetscInt c, c2;
5645733454dSMatthew G. Knepley 
5655733454dSMatthew G. Knepley     for (c = 0; c < 5; ++c)
5665733454dSMatthew G. Knepley       if (!used[c]) break;
5675733454dSMatthew G. Knepley     PetscCheck(c < 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell);
5685733454dSMatthew G. Knepley     // Match first edge to 3rd edge in newCone[2]
5695733454dSMatthew G. Knepley     {
5705733454dSMatthew G. Knepley       const PetscInt *fcone2, *fornt2, *farr2;
5715733454dSMatthew G. Knepley 
5725733454dSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt));
5735733454dSMatthew G. Knepley       farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]);
5745733454dSMatthew G. Knepley       PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2));
5755733454dSMatthew G. Knepley       farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]);
5765733454dSMatthew G. Knepley 
5775733454dSMatthew G. Knepley       const PetscInt e    = 2;
5785733454dSMatthew G. Knepley       const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]);
5795733454dSMatthew G. Knepley       // Trying to match edge `edge` with final orientation `eornt` of face `newCone[2]` to some edge of face `cone[c]` with orientation `ornt[c]`
5805733454dSMatthew G. Knepley       for (c2 = 0; c2 < 3; ++c2) {
5815733454dSMatthew G. Knepley         const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]);
5825733454dSMatthew G. Knepley         // Trying to match edge `edge2` with final orientation `eornt2`
5835733454dSMatthew G. Knepley         if (edge == edge2) {
5845733454dSMatthew G. Knepley           PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge);
5855733454dSMatthew G. Knepley           // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge`
5865733454dSMatthew G. Knepley           break;
5875733454dSMatthew G. Knepley         }
5885733454dSMatthew G. Knepley       }
5895733454dSMatthew G. Knepley       PetscCheck(c2 < 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in");
5905733454dSMatthew G. Knepley     }
5915733454dSMatthew G. Knepley     newCone[faces[4]] = cone[c];
5925733454dSMatthew G. Knepley     // Compute new orientation of face based on which edge was matched
5935733454dSMatthew G. Knepley     newOrnt[faces[4]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, c2, ornt[c]);
5945733454dSMatthew G. Knepley     PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt));
5955733454dSMatthew G. Knepley   }
5965733454dSMatthew G. Knepley   PetscCall(DMPlexSetCone(dm, cell, newCone));
5975733454dSMatthew G. Knepley   PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt));
5985733454dSMatthew G. Knepley   PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt));
5995733454dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
6005733454dSMatthew G. Knepley }
6015733454dSMatthew G. Knepley 
ReorderCell(PetscViewer viewer,DM dm,PetscInt cell,DMPolytopeType ct)602fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderCell(PetscViewer viewer, DM dm, PetscInt cell, DMPolytopeType ct)
603c7cc6ee6SMatthew G. Knepley {
604c7cc6ee6SMatthew G. Knepley   PetscFunctionBegin;
605c7cc6ee6SMatthew G. Knepley   switch (ct) {
606c7cc6ee6SMatthew G. Knepley   case DM_POLYTOPE_TRIANGLE:
607c7cc6ee6SMatthew G. Knepley   case DM_POLYTOPE_QUADRILATERAL:
608c7cc6ee6SMatthew G. Knepley     PetscCall(ReorderPolygon(dm, cell));
609c7cc6ee6SMatthew G. Knepley     break;
610c7cc6ee6SMatthew G. Knepley   case DM_POLYTOPE_TETRAHEDRON:
611fc8f8d65SMatthew G. Knepley     PetscCall(ReorderTetrahedron(viewer, dm, cell));
612c7cc6ee6SMatthew G. Knepley     break;
613c7cc6ee6SMatthew G. Knepley   case DM_POLYTOPE_HEXAHEDRON:
614c7cc6ee6SMatthew G. Knepley     PetscCall(ReorderHexahedron(dm, cell));
615c7cc6ee6SMatthew G. Knepley     break;
6165733454dSMatthew G. Knepley   case DM_POLYTOPE_TRI_PRISM:
6175733454dSMatthew G. Knepley     PetscCall(ReorderWedge(dm, cell));
6185733454dSMatthew G. Knepley     break;
619c7cc6ee6SMatthew G. Knepley   default:
620c7cc6ee6SMatthew G. Knepley     PetscCheck(0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Celltype %s is unsupported", DMPolytopeTypes[ct]);
621c7cc6ee6SMatthew G. Knepley     break;
622c7cc6ee6SMatthew G. Knepley   }
623c7cc6ee6SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
624c7cc6ee6SMatthew G. Knepley }
625c7cc6ee6SMatthew G. Knepley 
GetNumCellFaces(int nd,PetscInt * numCellFaces,DMPolytopeType * ct)6265733454dSMatthew G. Knepley static PetscErrorCode GetNumCellFaces(int nd, PetscInt *numCellFaces, DMPolytopeType *ct)
6275733454dSMatthew G. Knepley {
6285733454dSMatthew G. Knepley   PetscFunctionBegin;
6295733454dSMatthew G. Knepley   *ct = DM_POLYTOPE_POINT;
6305733454dSMatthew G. Knepley   switch (nd) {
6315733454dSMatthew G. Knepley   case 0:
6325733454dSMatthew G. Knepley     *numCellFaces = PETSC_DETERMINE;
6335733454dSMatthew G. Knepley     break;
6345733454dSMatthew G. Knepley   case 1:
6355733454dSMatthew G. Knepley     *numCellFaces = 3;
6365733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_TRIANGLE;
6375733454dSMatthew G. Knepley     break;
6385733454dSMatthew G. Knepley   case 2:
6395733454dSMatthew G. Knepley     *numCellFaces = 4;
6405733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_TETRAHEDRON;
6415733454dSMatthew G. Knepley     break;
6425733454dSMatthew G. Knepley   case 3:
6435733454dSMatthew G. Knepley     *numCellFaces = 4;
6445733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_QUADRILATERAL;
6455733454dSMatthew G. Knepley     break;
6465733454dSMatthew G. Knepley   case 4:
6475733454dSMatthew G. Knepley     *numCellFaces = 6;
6485733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_HEXAHEDRON;
6495733454dSMatthew G. Knepley     break;
6505733454dSMatthew G. Knepley   case 5:
6515733454dSMatthew G. Knepley     *numCellFaces = 5;
6525733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_PYRAMID;
6535733454dSMatthew G. Knepley     break;
6545733454dSMatthew G. Knepley   case 6:
6555733454dSMatthew G. Knepley     *numCellFaces = 5;
6565733454dSMatthew G. Knepley     *ct           = DM_POLYTOPE_TRI_PRISM;
6575733454dSMatthew G. Knepley     break;
6585733454dSMatthew G. Knepley   default:
6595733454dSMatthew G. Knepley     *numCellFaces = PETSC_DETERMINE;
6605733454dSMatthew G. Knepley   }
6615733454dSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
6625733454dSMatthew G. Knepley }
6635733454dSMatthew G. Knepley 
6642f0bd6dcSMichael Lange /*@C
6651d27aa22SBarry Smith   DMPlexCreateFluent - Create a `DMPLEX` mesh from a Fluent mesh file <http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm>.
6662f0bd6dcSMichael Lange 
667d083f849SBarry Smith   Collective
6682f0bd6dcSMichael Lange 
6692f0bd6dcSMichael Lange   Input Parameters:
6702f0bd6dcSMichael Lange + comm        - The MPI communicator
671a1cb98faSBarry Smith . viewer      - The `PetscViewer` associated with a Fluent mesh file
6722f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh
6732f0bd6dcSMichael Lange 
6742f0bd6dcSMichael Lange   Output Parameter:
675a1cb98faSBarry Smith . dm - The `DM` object representing the mesh
6762f0bd6dcSMichael Lange 
6772f0bd6dcSMichael Lange   Level: beginner
6782f0bd6dcSMichael Lange 
6791cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()`
6802f0bd6dcSMichael Lange @*/
DMPlexCreateFluent(MPI_Comm comm,PetscViewer viewer,PetscBool interpolate,DM * dm)681d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
682d71ae5a4SJacob Faibussowitsch {
683c7cc6ee6SMatthew G. Knepley   PetscInt        dim          = PETSC_DETERMINE;
684c7cc6ee6SMatthew G. Knepley   PetscInt        numCells     = 0;
685c7cc6ee6SMatthew G. Knepley   PetscInt        numVertices  = 0;
6865733454dSMatthew G. Knepley   PetscInt       *cellSizes    = NULL;
6875733454dSMatthew G. Knepley   DMPolytopeType *cellTypes    = NULL;
688c7cc6ee6SMatthew G. Knepley   PetscInt        numFaces     = 0;
689c7cc6ee6SMatthew G. Knepley   PetscInt       *faces        = NULL;
6905733454dSMatthew G. Knepley   PetscInt       *faceSizes    = NULL;
6915733454dSMatthew G. Knepley   PetscInt       *faceAdjCell  = NULL;
692c7cc6ee6SMatthew G. Knepley   PetscInt       *cellVertices = NULL;
693c7cc6ee6SMatthew G. Knepley   unsigned int   *faceZoneIDs  = NULL;
6947368db69SLisandro Dalcin   DMLabel         faceSets     = NULL;
695c7cc6ee6SMatthew G. Knepley   DMLabel        *zoneLabels   = NULL;
696c7cc6ee6SMatthew G. Knepley   const char    **zoneNames    = NULL;
697c7cc6ee6SMatthew G. Knepley   unsigned int    maxZoneID    = 0;
698c7cc6ee6SMatthew G. Knepley   PetscScalar    *coordsIn     = NULL;
699c7cc6ee6SMatthew G. Knepley   PetscScalar    *coords;
7003f6dc66eSMichael Lange   PetscSection    coordSection;
7013f6dc66eSMichael Lange   Vec             coordinates;
7025733454dSMatthew G. Knepley   PetscInt        coordSize, maxFaceSize = 0, totFaceVert = 0, f;
703c7cc6ee6SMatthew G. Knepley   PetscMPIInt     rank;
7042f0bd6dcSMichael Lange 
7052f0bd6dcSMichael Lange   PetscFunctionBegin;
7069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7072f0bd6dcSMichael Lange 
708dd400576SPatrick Sanan   if (rank == 0) {
7092f0bd6dcSMichael Lange     FluentSection s;
7105733454dSMatthew G. Knepley 
7115733454dSMatthew G. Knepley     s.data   = NULL;
712c7cc6ee6SMatthew G. Knepley     numFaces = PETSC_DETERMINE;
7132f0bd6dcSMichael Lange     do {
7149566063dSJacob Faibussowitsch       PetscCall(DMPlexCreateFluent_ReadSection(viewer, &s));
7152f0bd6dcSMichael Lange       if (s.index == 2) { /* Dimension */
716f7320561SMichael Lange         dim = s.nd;
717c7cc6ee6SMatthew G. Knepley         PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found dimension: %" PetscInt_FMT "\n", dim));
71819d58f9dSMichael Lange       } else if (s.index == 10 || s.index == 2010) { /* Vertices */
719c7cc6ee6SMatthew G. Knepley         if (s.zoneID == 0) {
720c7cc6ee6SMatthew G. Knepley           numVertices = s.last;
721c7cc6ee6SMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of vertices: %" PetscInt_FMT "\n", numVertices));
722c7cc6ee6SMatthew G. Knepley         } else {
723c7cc6ee6SMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found vertex coordinates\n"));
72428b400f6SJacob Faibussowitsch           PetscCheck(!coordsIn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files");
7251a9c30ecSMatthew G. Knepley           coordsIn = (PetscScalar *)s.data;
7263f6dc66eSMichael Lange         }
7272f0bd6dcSMichael Lange 
72819d58f9dSMichael Lange       } else if (s.index == 12 || s.index == 2012) { /* Cells */
729c7cc6ee6SMatthew G. Knepley         if (s.zoneID == 0) {
730c7cc6ee6SMatthew G. Knepley           numCells = s.last;
731c7cc6ee6SMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cells %" PetscInt_FMT "\n", numCells));
732c7cc6ee6SMatthew G. Knepley         } else {
7335733454dSMatthew G. Knepley           PetscCall(PetscMalloc2(numCells, &cellSizes, numCells, &cellTypes));
7349de2bd2cSJose E. Roman           for (PetscInt c = 0; c < numCells; ++c) PetscCall(GetNumCellFaces(s.nd ? s.nd : (int)((PetscInt *)s.data)[c], &cellSizes[c], &cellTypes[c]));
7355733454dSMatthew G. Knepley           PetscCall(PetscFree(s.data));
7365733454dSMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cell faces %" PetscInt_FMT "\n", numCells && s.nd ? cellSizes[0] : 0));
737f7320561SMichael Lange         }
73819d58f9dSMichael Lange       } else if (s.index == 13 || s.index == 2013) { /* Facets */
739f7320561SMichael Lange         if (s.zoneID == 0) {                         /* Header section */
740930bae4bSMatthew G. Knepley           numFaces = (PetscInt)(s.last - s.first + 1);
7415733454dSMatthew G. Knepley           PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of faces %" PetscInt_FMT " face vertices: %d\n", numFaces, s.nd));
742f7320561SMichael Lange         } else { /* Data section */
7435733454dSMatthew G. Knepley           PetscInt *tmp;
7445733454dSMatthew G. Knepley           PetscInt  totSize = 0, offset = 0, doffset;
745cf554e2cSMatthew G. Knepley 
74608401ef6SPierre Jolivet           PetscCheck(numFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file");
7475733454dSMatthew G. Knepley           if (!faceZoneIDs) PetscCall(PetscMalloc3(numFaces, &faceSizes, numFaces * 2, &faceAdjCell, numFaces, &faceZoneIDs));
7485733454dSMatthew G. Knepley           // Record the zoneID and face size for each face set
7495733454dSMatthew G. Knepley           for (unsigned int z = s.first - 1; z < s.last; z++) {
7505733454dSMatthew G. Knepley             faceZoneIDs[z] = s.zoneID;
7515733454dSMatthew G. Knepley             if (s.nd) {
7525733454dSMatthew G. Knepley               faceSizes[z] = s.nd;
7535733454dSMatthew G. Knepley             } else {
7545733454dSMatthew G. Knepley               faceSizes[z] = ((PetscInt *)s.data)[offset];
7555733454dSMatthew G. Knepley               offset += faceSizes[z] + 3;
7565733454dSMatthew G. Knepley             }
7575733454dSMatthew G. Knepley             totSize += faceSizes[z];
7585733454dSMatthew G. Knepley             maxFaceSize = PetscMax(maxFaceSize, faceSizes[z]);
7595733454dSMatthew G. Knepley           }
7605733454dSMatthew G. Knepley 
7615733454dSMatthew G. Knepley           offset  = totFaceVert;
7625733454dSMatthew G. Knepley           doffset = s.nd ? 0 : 1;
7635733454dSMatthew G. Knepley           PetscCall(PetscMalloc1(totFaceVert + totSize, &tmp));
7645733454dSMatthew G. Knepley           if (faces) PetscCall(PetscArraycpy(tmp, faces, totFaceVert));
7655733454dSMatthew G. Knepley           PetscCall(PetscFree(faces));
7665733454dSMatthew G. Knepley           totFaceVert += totSize;
7675733454dSMatthew G. Knepley           faces = tmp;
7685733454dSMatthew G. Knepley 
7695733454dSMatthew G. Knepley           // Record face vertices and adjacent faces
7705733454dSMatthew G. Knepley           const PetscInt Nfz = s.last - s.first + 1;
7715733454dSMatthew G. Knepley           for (PetscInt f = 0; f < Nfz; ++f) {
7725733454dSMatthew G. Knepley             const PetscInt face     = f + s.first - 1;
7735733454dSMatthew G. Knepley             const PetscInt faceSize = faceSizes[face];
7745733454dSMatthew G. Knepley 
7755733454dSMatthew G. Knepley             for (PetscInt v = 0; v < faceSize; ++v) faces[offset + v] = ((PetscInt *)s.data)[doffset + v];
7765733454dSMatthew G. Knepley             faceAdjCell[face * 2 + 0] = ((PetscInt *)s.data)[doffset + faceSize + 0];
7775733454dSMatthew G. Knepley             faceAdjCell[face * 2 + 1] = ((PetscInt *)s.data)[doffset + faceSize + 1];
7785733454dSMatthew G. Knepley             offset += faceSize;
7795733454dSMatthew G. Knepley             doffset += faceSize + (s.nd ? 2 : 3);
7805733454dSMatthew G. Knepley           }
7819566063dSJacob Faibussowitsch           PetscCall(PetscFree(s.data));
782f7320561SMichael Lange         }
783c7cc6ee6SMatthew G. Knepley       } else if (s.index == 39) { /* Label information */
784c7cc6ee6SMatthew G. Knepley         if (s.zoneID >= maxZoneID) {
785c7cc6ee6SMatthew G. Knepley           DMLabel     *tmpL;
786c7cc6ee6SMatthew G. Knepley           const char **tmp;
787c7cc6ee6SMatthew G. Knepley           unsigned int newmax = maxZoneID + 1;
788c7cc6ee6SMatthew G. Knepley 
789c7cc6ee6SMatthew G. Knepley           while (newmax < s.zoneID + 1) newmax *= 2;
790c7cc6ee6SMatthew G. Knepley           PetscCall(PetscCalloc2(newmax, &tmp, newmax, &tmpL));
791c7cc6ee6SMatthew G. Knepley           for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) {
792c7cc6ee6SMatthew G. Knepley             tmp[i]  = zoneNames[i];
793c7cc6ee6SMatthew G. Knepley             tmpL[i] = zoneLabels[i];
794c7cc6ee6SMatthew G. Knepley           }
795c7cc6ee6SMatthew G. Knepley           maxZoneID = newmax;
796c7cc6ee6SMatthew G. Knepley           PetscCall(PetscFree2(zoneNames, zoneLabels));
797c7cc6ee6SMatthew G. Knepley           zoneNames  = tmp;
798c7cc6ee6SMatthew G. Knepley           zoneLabels = tmpL;
799c7cc6ee6SMatthew G. Knepley         }
800c7cc6ee6SMatthew G. Knepley         zoneNames[s.zoneID] = (const char *)s.data;
8012f0bd6dcSMichael Lange       }
8022f0bd6dcSMichael Lange     } while (s.index >= 0);
8032f0bd6dcSMichael Lange   }
8049566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm));
80508401ef6SPierre Jolivet   PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension");
806f7320561SMichael Lange 
807f7320561SMichael Lange   /* Allocate cell-vertex mesh */
8089566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
8099566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
8109566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(*dm, dim));
8115733454dSMatthew G. Knepley   // We do not want this label automatically computed, instead we fill it here
8125733454dSMatthew G. Knepley   PetscCall(DMCreateLabel(*dm, "celltype"));
813c7cc6ee6SMatthew G. Knepley   PetscCall(DMPlexSetChart(*dm, 0, numCells + numFaces + numVertices));
814dd400576SPatrick Sanan   if (rank == 0) {
81508401ef6SPierre Jolivet     PetscCheck(numCells >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file");
8165733454dSMatthew G. Knepley     for (PetscInt c = 0; c < numCells; ++c) {
8175733454dSMatthew G. Knepley       PetscCall(DMPlexSetConeSize(*dm, c, cellSizes[c]));
8185733454dSMatthew G. Knepley       PetscCall(DMPlexSetCellType(*dm, c, cellTypes[c]));
819c7cc6ee6SMatthew G. Knepley     }
8205733454dSMatthew G. Knepley     for (PetscInt v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
8215733454dSMatthew G. Knepley     for (PetscInt f = 0; f < numFaces; ++f) {
8225733454dSMatthew G. Knepley       DMPolytopeType ct;
8235733454dSMatthew G. Knepley 
8245733454dSMatthew G. Knepley       switch (faceSizes[f]) {
8255733454dSMatthew G. Knepley       case 2:
8265733454dSMatthew G. Knepley         ct = DM_POLYTOPE_SEGMENT;
8275733454dSMatthew G. Knepley         break;
8285733454dSMatthew G. Knepley       case 3:
8295733454dSMatthew G. Knepley         ct = DM_POLYTOPE_TRIANGLE;
8305733454dSMatthew G. Knepley         break;
8315733454dSMatthew G. Knepley       case 4:
8325733454dSMatthew G. Knepley         ct = DM_POLYTOPE_QUADRILATERAL;
8335733454dSMatthew G. Knepley         break;
8345733454dSMatthew G. Knepley       default:
8355733454dSMatthew G. Knepley         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file with cone size %" PetscInt_FMT, faceSizes[f]);
8365733454dSMatthew G. Knepley       }
8375733454dSMatthew G. Knepley       PetscCall(DMPlexSetConeSize(*dm, f + numCells + numVertices, faceSizes[f]));
8385733454dSMatthew G. Knepley       PetscCall(DMPlexSetCellType(*dm, f + numCells + numVertices, ct));
8395733454dSMatthew G. Knepley     }
840f7320561SMichael Lange   }
8419566063dSJacob Faibussowitsch   PetscCall(DMSetUp(*dm));
842f7320561SMichael Lange 
843dd400576SPatrick Sanan   if (rank == 0 && faces) {
8445733454dSMatthew G. Knepley     PetscSection s;
8455733454dSMatthew G. Knepley     PetscInt    *cones, csize, foffset = 0;
846f7320561SMichael Lange 
847c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexGetCones(*dm, &cones));
8485733454dSMatthew G. Knepley     PetscCall(DMPlexGetConeSection(*dm, &s));
8495733454dSMatthew G. Knepley     PetscCall(PetscSectionGetConstrainedStorageSize(s, &csize));
8505733454dSMatthew G. Knepley     for (PetscInt c = 0; c < csize; ++c) cones[c] = -1;
851c7cc6ee6SMatthew G. Knepley     for (PetscInt f = 0; f < numFaces; f++) {
8525733454dSMatthew G. Knepley       const PetscInt cl   = faceAdjCell[f * 2 + 0] - 1;
8535733454dSMatthew G. Knepley       const PetscInt cr   = faceAdjCell[f * 2 + 1] - 1;
854c7cc6ee6SMatthew G. Knepley       const PetscInt face = f + numCells + numVertices;
855c7cc6ee6SMatthew G. Knepley       PetscInt       fcone[16];
856c7cc6ee6SMatthew G. Knepley 
857c7cc6ee6SMatthew G. Knepley       // How could Fluent define the outward normal differently? Is there no end to the pain?
858c7cc6ee6SMatthew G. Knepley       if (dim == 3) {
859c7cc6ee6SMatthew G. Knepley         if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, -1));
860c7cc6ee6SMatthew G. Knepley         if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, 0));
861c7cc6ee6SMatthew G. Knepley       } else {
862c7cc6ee6SMatthew G. Knepley         if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, 0));
863c7cc6ee6SMatthew G. Knepley         if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, -1));
8649371c9d4SSatish Balay       }
8655733454dSMatthew G. Knepley       PetscCheck(faceSizes[f] < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of face vertices %" PetscInt_FMT " exceeds temporary storage", faceSizes[f]);
8665733454dSMatthew G. Knepley       for (PetscInt v = 0; v < faceSizes[f]; ++v) fcone[v] = faces[foffset + v] + numCells - 1;
8675733454dSMatthew G. Knepley       foffset += faceSizes[f];
868c7cc6ee6SMatthew G. Knepley       PetscCall(DMPlexSetCone(*dm, face, fcone));
869f7320561SMichael Lange     }
870f7320561SMichael Lange   }
8719566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
8729566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
873c7cc6ee6SMatthew G. Knepley   if (dim == 3) {
8745fd9971aSMatthew G. Knepley     DM idm;
875f7320561SMichael Lange 
876c7cc6ee6SMatthew G. Knepley     PetscCall(DMCreate(PetscObjectComm((PetscObject)*dm), &idm));
877c7cc6ee6SMatthew G. Knepley     PetscCall(DMSetType(idm, DMPLEX));
878c7cc6ee6SMatthew G. Knepley     PetscCall(DMSetDimension(idm, dim));
879c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexInterpolateFaces_Internal(*dm, 1, idm));
8809566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
881f7320561SMichael Lange     *dm = idm;
882f7320561SMichael Lange   }
883c7cc6ee6SMatthew G. Knepley   PetscCall(DMViewFromOptions(*dm, NULL, "-cas_dm_view"));
884c7cc6ee6SMatthew G. Knepley   if (rank == 0 && faces) {
8855733454dSMatthew G. Knepley     for (PetscInt c = 0; c < numCells; ++c) PetscCall(ReorderCell(viewer, *dm, c, cellTypes[c]));
886c7cc6ee6SMatthew G. Knepley   }
887f7320561SMichael Lange 
888dd400576SPatrick Sanan   if (rank == 0 && faces) {
8895733454dSMatthew G. Knepley     PetscInt        joinSize, meetSize, *fverts, cells[2];
890631eb916SMichael Lange     const PetscInt *join, *meet;
8915733454dSMatthew G. Knepley     PetscInt        foffset = 0;
8925733454dSMatthew G. Knepley 
8935733454dSMatthew G. Knepley     PetscCall(PetscMalloc1(maxFaceSize, &fverts));
894ec78a56aSMichael Lange     /* Mark facets by finding the full join of all adjacent vertices */
895ec78a56aSMichael Lange     for (f = 0; f < numFaces; f++) {
8965733454dSMatthew G. Knepley       const PetscInt cl = faceAdjCell[f * 2 + 0] - 1;
8975733454dSMatthew G. Knepley       const PetscInt cr = faceAdjCell[f * 2 + 1] - 1;
898c7cc6ee6SMatthew G. Knepley       const PetscInt id = (PetscInt)faceZoneIDs[f];
899c7cc6ee6SMatthew G. Knepley 
900631eb916SMichael Lange       if (cl > 0 && cr > 0) {
901631eb916SMichael Lange         /* If we know both adjoining cells we can use a single-level meet */
9029371c9d4SSatish Balay         cells[0] = cl;
9039371c9d4SSatish Balay         cells[1] = cr;
9049566063dSJacob Faibussowitsch         PetscCall(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet));
905c7cc6ee6SMatthew G. Knepley         PetscCheck(meetSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT " cells: %" PetscInt_FMT ", %" PetscInt_FMT, f, cl, cr);
906c7cc6ee6SMatthew G. Knepley         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], id));
907c7cc6ee6SMatthew G. Knepley         if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], meet[0], 1));
9085733454dSMatthew G. Knepley         PetscCall(DMPlexRestoreMeet(*dm, meetSize, fverts, &meetSize, &meet));
909631eb916SMichael Lange       } else {
9105733454dSMatthew G. Knepley         for (PetscInt fi = 0; fi < faceSizes[f]; fi++) fverts[fi] = faces[foffset + fi] + numCells - 1;
9115733454dSMatthew G. Knepley         PetscCall(DMPlexGetFullJoin(*dm, faceSizes[f], fverts, &joinSize, &join));
91263a3b9bcSJacob Faibussowitsch         PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f);
913c7cc6ee6SMatthew G. Knepley         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], id));
914c7cc6ee6SMatthew G. Knepley         if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], join[0], 1));
9155733454dSMatthew G. Knepley         PetscCall(DMPlexRestoreJoin(*dm, joinSize, fverts, &joinSize, &join));
916ec78a56aSMichael Lange       }
9175733454dSMatthew G. Knepley       foffset += faceSizes[f];
918631eb916SMichael Lange     }
9199566063dSJacob Faibussowitsch     PetscCall(PetscFree(fverts));
920ec78a56aSMichael Lange   }
921ec78a56aSMichael Lange 
9227368db69SLisandro Dalcin   { /* Create Face Sets label at all processes */
9239371c9d4SSatish Balay     enum {
9249371c9d4SSatish Balay       n = 1
9259371c9d4SSatish Balay     };
9267368db69SLisandro Dalcin     PetscBool flag[n];
9277368db69SLisandro Dalcin 
9287368db69SLisandro Dalcin     flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE;
929*5440e5dcSBarry Smith     PetscCallMPI(MPI_Bcast(flag, n, MPI_C_BOOL, 0, comm));
9309566063dSJacob Faibussowitsch     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
931c7cc6ee6SMatthew G. Knepley     // TODO Code to create all the zone labels on each process
932c7cc6ee6SMatthew G. Knepley   }
933c7cc6ee6SMatthew G. Knepley 
934c7cc6ee6SMatthew G. Knepley   if (!interpolate) {
935c7cc6ee6SMatthew G. Knepley     DM udm;
936c7cc6ee6SMatthew G. Knepley 
937c7cc6ee6SMatthew G. Knepley     PetscCall(DMPlexUninterpolate(*dm, &udm));
938c7cc6ee6SMatthew G. Knepley     PetscCall(DMDestroy(dm));
939c7cc6ee6SMatthew G. Knepley     *dm = udm;
9407368db69SLisandro Dalcin   }
9417368db69SLisandro Dalcin 
9423f6dc66eSMichael Lange   /* Read coordinates */
9439566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
9449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
9459566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim));
9469566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
947c7cc6ee6SMatthew G. Knepley   for (PetscInt v = numCells; v < numCells + numVertices; ++v) {
9489566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, dim));
9499566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim));
9503f6dc66eSMichael Lange   }
9519566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
9529566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
9539566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
9549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
9559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
9569566063dSJacob Faibussowitsch   PetscCall(VecSetType(coordinates, VECSTANDARD));
9579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
958dd400576SPatrick Sanan   if (rank == 0 && coordsIn) {
959c7cc6ee6SMatthew G. Knepley     for (PetscInt v = 0; v < numVertices; ++v) {
960c7cc6ee6SMatthew G. Knepley       for (PetscInt d = 0; d < dim; ++d) coords[v * dim + d] = coordsIn[v * dim + d];
9613f6dc66eSMichael Lange     }
9623f6dc66eSMichael Lange   }
9639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
9649566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
9659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
9667368db69SLisandro Dalcin 
967dd400576SPatrick Sanan   if (rank == 0) {
9689566063dSJacob Faibussowitsch     PetscCall(PetscFree(cellVertices));
9695733454dSMatthew G. Knepley     PetscCall(PetscFree2(cellSizes, cellTypes));
9709566063dSJacob Faibussowitsch     PetscCall(PetscFree(faces));
9715733454dSMatthew G. Knepley     PetscCall(PetscFree3(faceSizes, faceAdjCell, faceZoneIDs));
9729566063dSJacob Faibussowitsch     PetscCall(PetscFree(coordsIn));
973c7cc6ee6SMatthew G. Knepley     if (zoneNames)
974c7cc6ee6SMatthew G. Knepley       for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) PetscCall(PetscFree(zoneNames[i]));
975c7cc6ee6SMatthew G. Knepley     PetscCall(PetscFree2(zoneNames, zoneLabels));
976f7320561SMichael Lange   }
9773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9782f0bd6dcSMichael Lange }
979