xref: /petsc/src/dm/impls/plex/plexfluent.c (revision 5c0a1baa0f4a9e6dc65d5aad147407b3395fcdb1)
12f0bd6dcSMichael Lange #define PETSCDM_DLL
22f0bd6dcSMichael Lange #include <petsc-private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
32f0bd6dcSMichael Lange 
42f0bd6dcSMichael Lange #undef __FUNCT__
52f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluentFromFile"
62f0bd6dcSMichael Lange /*@C
72f0bd6dcSMichael Lange   DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file
82f0bd6dcSMichael Lange 
92f0bd6dcSMichael Lange + comm        - The MPI communicator
102f0bd6dcSMichael Lange . filename    - Name of the Fluent mesh file
112f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh
122f0bd6dcSMichael Lange 
132f0bd6dcSMichael Lange   Output Parameter:
142f0bd6dcSMichael Lange . dm  - The DM object representing the mesh
152f0bd6dcSMichael Lange 
162f0bd6dcSMichael Lange   Level: beginner
172f0bd6dcSMichael Lange 
182f0bd6dcSMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate()
192f0bd6dcSMichael Lange @*/
202f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
212f0bd6dcSMichael Lange {
222f0bd6dcSMichael Lange   PetscViewer     viewer;
232f0bd6dcSMichael Lange   PetscErrorCode  ierr;
242f0bd6dcSMichael Lange 
252f0bd6dcSMichael Lange   PetscFunctionBegin;
262f0bd6dcSMichael Lange   /* Create file viewer and build plex */
272f0bd6dcSMichael Lange   ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr);
282f0bd6dcSMichael Lange   ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr);
292f0bd6dcSMichael Lange   ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr);
302f0bd6dcSMichael Lange   ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr);
312f0bd6dcSMichael Lange   ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr);
322f0bd6dcSMichael Lange   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
332f0bd6dcSMichael Lange   PetscFunctionReturn(0);
342f0bd6dcSMichael Lange }
352f0bd6dcSMichael Lange 
362f0bd6dcSMichael Lange #undef __FUNCT__
372f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent_ReadString"
382f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim)
392f0bd6dcSMichael Lange {
402f0bd6dcSMichael Lange   PetscInt i = 0;
412f0bd6dcSMichael Lange   PetscErrorCode ierr;
422f0bd6dcSMichael Lange 
432f0bd6dcSMichael Lange   PetscFunctionBegin;
442f0bd6dcSMichael Lange   do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, PETSC_CHAR);CHKERRQ(ierr);}
452f0bd6dcSMichael Lange   while (buffer[i-1] != '\0' && buffer[i-1] != delim);
462f0bd6dcSMichael Lange   buffer[i] = '\0';
472f0bd6dcSMichael Lange   PetscFunctionReturn(0);
482f0bd6dcSMichael Lange }
492f0bd6dcSMichael Lange 
502f0bd6dcSMichael Lange #undef __FUNCT__
5119d58f9dSMichael Lange #define __FUNCT__ "DMPlexCreateFluent_ReadValues"
5219d58f9dSMichael Lange PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary)
5319d58f9dSMichael Lange {
5419d58f9dSMichael Lange   int            fdes;
5519d58f9dSMichael Lange   FILE          *file;
5619d58f9dSMichael Lange   PetscInt       i;
5719d58f9dSMichael Lange   PetscErrorCode ierr;
5819d58f9dSMichael Lange 
5919d58f9dSMichael Lange   PetscFunctionBegin;
6019d58f9dSMichael Lange   if (binary) {
6119d58f9dSMichael Lange     /* Extract raw file descriptor to read binary block */
6219d58f9dSMichael Lange     ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr);
6319d58f9dSMichael Lange     fflush(file); fdes = fileno(file);
6419d58f9dSMichael Lange   }
6519d58f9dSMichael Lange 
6619d58f9dSMichael Lange   if (!binary && dtype == PETSC_INT) {
6719d58f9dSMichael Lange     char cbuf[256];
6819d58f9dSMichael Lange     int ibuf, snum;
6919d58f9dSMichael Lange     /* Parse hexadecimal ascii integers */
7019d58f9dSMichael Lange     for (i = 0; i < count; i++) {
7119d58f9dSMichael Lange       ierr = PetscViewerRead(viewer, cbuf, 1, PETSC_STRING);CHKERRQ(ierr);
7219d58f9dSMichael Lange       snum = sscanf(cbuf, "%x", &ibuf);
7319d58f9dSMichael Lange       if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
7419d58f9dSMichael Lange       ((PetscInt*)data)[i] = (PetscInt)ibuf;
7519d58f9dSMichael Lange     }
7619d58f9dSMichael Lange   } else if (binary && dtype == PETSC_INT) {
7719d58f9dSMichael Lange     /* Always read 32-bit ints and cast to PetscInt */
7819d58f9dSMichael Lange     int *ibuf;
7919d58f9dSMichael Lange     ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr);
8019d58f9dSMichael Lange     ierr = PetscBinaryRead(fdes, ibuf, count, PETSC_ENUM);CHKERRQ(ierr);
8119d58f9dSMichael Lange     ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr);
8219d58f9dSMichael Lange     for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]);
8319d58f9dSMichael Lange     ierr = PetscFree(ibuf);CHKERRQ(ierr);
8419d58f9dSMichael Lange 
8519d58f9dSMichael Lange  } else if (binary && dtype == PETSC_SCALAR) {
8619d58f9dSMichael Lange     float *fbuf;
8719d58f9dSMichael Lange     /* Always read 32-bit floats and cast to PetscScalar */
8819d58f9dSMichael Lange     ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr);
8919d58f9dSMichael Lange     ierr = PetscBinaryRead(fdes, fbuf, count, PETSC_FLOAT);CHKERRQ(ierr);
9019d58f9dSMichael Lange     ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr);
9119d58f9dSMichael Lange     for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]);
9219d58f9dSMichael Lange     ierr = PetscFree(fbuf);CHKERRQ(ierr);
9319d58f9dSMichael Lange   } else {
9419d58f9dSMichael Lange     ierr = PetscViewerASCIIRead(viewer, data, count, dtype);CHKERRQ(ierr);
9519d58f9dSMichael Lange   }
9619d58f9dSMichael Lange   PetscFunctionReturn(0);
9719d58f9dSMichael Lange }
9819d58f9dSMichael Lange 
9919d58f9dSMichael Lange #undef __FUNCT__
1002f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent_ReadSection"
1012f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s)
1022f0bd6dcSMichael Lange {
1032f0bd6dcSMichael Lange   char           buffer[PETSC_MAX_PATH_LEN];
1042f0bd6dcSMichael Lange   int            snum;
1052f0bd6dcSMichael Lange   PetscErrorCode ierr;
1062f0bd6dcSMichael Lange 
1072f0bd6dcSMichael Lange   PetscFunctionBegin;
1082f0bd6dcSMichael Lange   /* Fast-forward to next section and derive its index */
1092f0bd6dcSMichael Lange   ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
1102f0bd6dcSMichael Lange   ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr);
1112f0bd6dcSMichael Lange   snum = sscanf(buffer, "%d", &(s->index));
1122f0bd6dcSMichael Lange   /* If we can't match an index return -1 to signal end-of-file */
1132f0bd6dcSMichael Lange   if (snum < 1) {s->index = -1;   PetscFunctionReturn(0);}
1142f0bd6dcSMichael Lange 
1152f0bd6dcSMichael Lange   if (s->index == 0) {           /* Comment */
1162f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1172f0bd6dcSMichael Lange 
1182f0bd6dcSMichael Lange   } else if (s->index == 2) {    /* Dimension */
1192f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1202f0bd6dcSMichael Lange     snum = sscanf(buffer, "%d", &(s->nd));
1212f0bd6dcSMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1222f0bd6dcSMichael Lange 
12319d58f9dSMichael Lange   } else if (s->index == 10 || s->index == 2010) {   /* Vertices */
1242f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1252f0bd6dcSMichael Lange     snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
1262f0bd6dcSMichael Lange     if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1273f6dc66eSMichael Lange     if (s->zoneID > 0) {
1283f6dc66eSMichael Lange       PetscInt numCoords = s->last - s->first + 1;
1293f6dc66eSMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
1303f6dc66eSMichael Lange       ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr);
131*5c0a1baaSMatthew G. Knepley       ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
13219d58f9dSMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1333f6dc66eSMichael Lange     }
13419d58f9dSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1352f0bd6dcSMichael Lange 
13619d58f9dSMichael Lange   } else if (s->index == 12 || s->index == 2012) {   /* Cells */
1372f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1382f0bd6dcSMichael Lange     snum = sscanf(buffer, "(%x", &(s->zoneID));
1392f0bd6dcSMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1402f0bd6dcSMichael Lange     if (s->zoneID == 0) {  /* Header section */
1412f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
1422f0bd6dcSMichael Lange       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
143f7320561SMichael Lange     } else {               /* Data section */
1442f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
1452f0bd6dcSMichael Lange       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
146895fcc3eSMichael Lange       if (s->nd == 0) {
147895fcc3eSMichael Lange         /* Read cell type definitions for mixed cells */
14819d58f9dSMichael Lange         PetscInt numCells = s->last - s->first + 1;
149895fcc3eSMichael Lange         ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
15019d58f9dSMichael Lange         ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr);
151*5c0a1baaSMatthew G. Knepley         ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
15219d58f9dSMichael Lange         ierr = PetscFree(s->data);CHKERRQ(ierr);
153895fcc3eSMichael Lange         ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
154895fcc3eSMichael Lange       }
1552f0bd6dcSMichael Lange     }
1562f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1572f0bd6dcSMichael Lange 
15819d58f9dSMichael Lange   } else if (s->index == 13 || s->index == 2013) {   /* Faces */
1592f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1602f0bd6dcSMichael Lange     snum = sscanf(buffer, "(%x", &(s->zoneID));
1612f0bd6dcSMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
1622f0bd6dcSMichael Lange     if (s->zoneID == 0) {  /* Header section */
1632f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
1642f0bd6dcSMichael Lange       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
165f7320561SMichael Lange     } else {               /* Data section */
16619d58f9dSMichael Lange       PetscInt f, numEntries, numFaces, numFaceVert;
1672f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
1682f0bd6dcSMichael Lange       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
169f7320561SMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
170f7320561SMichael Lange       switch (s->nd) {
171895fcc3eSMichael Lange       case 0: numEntries = PETSC_DETERMINE; break;
172f7320561SMichael Lange       case 2: numEntries = 2 + 2; break;  /* linear */
173f7320561SMichael Lange       case 3: numEntries = 2 + 3; break;  /* triangular */
174f7320561SMichael Lange       case 4: numEntries = 2 + 4; break;  /* quadrilateral */
175f7320561SMichael Lange       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file");
176f7320561SMichael Lange       }
177f7320561SMichael Lange       numFaces = s->last-s->first + 1;
178895fcc3eSMichael Lange       if (numEntries != PETSC_DETERMINE) {
179895fcc3eSMichael Lange         /* Allocate space only if we already know the size of the block */
180f7320561SMichael Lange         ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr);
181895fcc3eSMichael Lange       }
182f7320561SMichael Lange       for (f = 0; f < numFaces; f++) {
183895fcc3eSMichael Lange         if (s->nd == 0) {
184895fcc3eSMichael Lange           /* Determine the size of the block for "mixed" facets */
185*5c0a1baaSMatthew G. Knepley           ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
186895fcc3eSMichael Lange           if (numEntries == PETSC_DETERMINE) {
187895fcc3eSMichael Lange             numEntries = numFaceVert + 2;
188895fcc3eSMichael Lange             ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr);
189895fcc3eSMichael Lange           } else {
190895fcc3eSMichael Lange             if (numEntries != numFaceVert + 2) {
191895fcc3eSMichael Lange               SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files");
192895fcc3eSMichael Lange             }
193895fcc3eSMichael Lange           }
194895fcc3eSMichael Lange         }
195*5c0a1baaSMatthew G. Knepley         ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr);
196f7320561SMichael Lange       }
197895fcc3eSMichael Lange       s->nd = numEntries - 2;
198f7320561SMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1992f0bd6dcSMichael Lange     }
2002f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
2012f0bd6dcSMichael Lange 
2022f0bd6dcSMichael Lange   } else {                       /* Unknown section type */
2032f0bd6dcSMichael Lange     PetscInt depth = 1;
2042f0bd6dcSMichael Lange     do {
2052f0bd6dcSMichael Lange       /* Match parentheses when parsing unknown sections */
2062f0bd6dcSMichael Lange       do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, PETSC_CHAR);CHKERRQ(ierr);}
2072f0bd6dcSMichael Lange       while (buffer[0] != '(' && buffer[0] != ')');
2082f0bd6dcSMichael Lange       if (buffer[0] == '(') depth++;
2092f0bd6dcSMichael Lange       if (buffer[0] == ')') depth--;
2102f0bd6dcSMichael Lange     } while (depth > 0);
2112f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr);
2122f0bd6dcSMichael Lange   }
2132f0bd6dcSMichael Lange   PetscFunctionReturn(0);
2142f0bd6dcSMichael Lange }
2152f0bd6dcSMichael Lange 
2162f0bd6dcSMichael Lange #undef __FUNCT__
2172f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent"
2182f0bd6dcSMichael Lange /*@C
2192f0bd6dcSMichael Lange   DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file.
2202f0bd6dcSMichael Lange 
2212f0bd6dcSMichael Lange   Collective on comm
2222f0bd6dcSMichael Lange 
2232f0bd6dcSMichael Lange   Input Parameters:
2242f0bd6dcSMichael Lange + comm  - The MPI communicator
2252f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file
2262f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh
2272f0bd6dcSMichael Lange 
2282f0bd6dcSMichael Lange   Output Parameter:
2292f0bd6dcSMichael Lange . dm  - The DM object representing the mesh
2302f0bd6dcSMichael Lange 
2312f0bd6dcSMichael Lange   Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm
2322f0bd6dcSMichael Lange 
2332f0bd6dcSMichael Lange   Level: beginner
2342f0bd6dcSMichael Lange 
2352f0bd6dcSMichael Lange .keywords: mesh, fluent, case
2362f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate()
2372f0bd6dcSMichael Lange @*/
2382f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
2392f0bd6dcSMichael Lange {
2402f0bd6dcSMichael Lange   PetscMPIInt    rank;
241f7320561SMichael Lange   PetscInt       c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE;
242f7320561SMichael Lange   PetscInt       numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE;
243ec78a56aSMichael Lange   PetscInt      *faces = NULL, *cellVertices, *faceZoneIDs = NULL;
2443f6dc66eSMichael Lange   PetscInt       d, coordSize;
2453f6dc66eSMichael Lange   PetscScalar   *coords, *coordsIn = NULL;
2463f6dc66eSMichael Lange   PetscSection   coordSection;
2473f6dc66eSMichael Lange   Vec            coordinates;
2482f0bd6dcSMichael Lange   PetscErrorCode ierr;
2492f0bd6dcSMichael Lange 
2502f0bd6dcSMichael Lange   PetscFunctionBegin;
2512f0bd6dcSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
2522f0bd6dcSMichael Lange 
2532f0bd6dcSMichael Lange   if (!rank) {
2542f0bd6dcSMichael Lange     FluentSection s;
2552f0bd6dcSMichael Lange     do {
2562f0bd6dcSMichael Lange       ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr);
2572f0bd6dcSMichael Lange       if (s.index == 2) {                 /* Dimension */
258f7320561SMichael Lange         dim = s.nd;
2592f0bd6dcSMichael Lange 
26019d58f9dSMichael Lange       } else if (s.index == 10 || s.index == 2010) { /* Vertices */
261f7320561SMichael Lange         if (s.zoneID == 0) numVertices = s.last;
2623f6dc66eSMichael Lange         else {
2633f6dc66eSMichael Lange           if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files");
2643f6dc66eSMichael Lange           coordsIn = s.data;
2653f6dc66eSMichael Lange         }
2662f0bd6dcSMichael Lange 
26719d58f9dSMichael Lange       } else if (s.index == 12 || s.index == 2012) { /* Cells */
268f7320561SMichael Lange         if (s.zoneID == 0) numCells = s.last;
269f7320561SMichael Lange         else {
270f7320561SMichael Lange           switch (s.nd) {
271895fcc3eSMichael Lange           case 0: numCellVertices = PETSC_DETERMINE; break;
272f7320561SMichael Lange           case 1: numCellVertices = 3; break;  /* triangular */
273f7320561SMichael Lange           case 2: numCellVertices = 4; break;  /* tetrahedral */
274f7320561SMichael Lange           case 3: numCellVertices = 4; break;  /* quadrilateral */
275f7320561SMichael Lange           case 4: numCellVertices = 8; break;  /* hexahedral */
276f7320561SMichael Lange           case 5: numCellVertices = 5; break;  /* pyramid */
277f7320561SMichael Lange           case 6: numCellVertices = 6; break;  /* wedge */
278895fcc3eSMichael Lange           default: numCellVertices = PETSC_DETERMINE;
279f7320561SMichael Lange           }
280f7320561SMichael Lange         }
2812f0bd6dcSMichael Lange 
28219d58f9dSMichael Lange       } else if (s.index == 13 || s.index == 2013) { /* Facets */
283f7320561SMichael Lange         if (s.zoneID == 0) {  /* Header section */
284f7320561SMichael Lange           numFaces = s.last - s.first + 1;
28535462f7fSMichael Lange           if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE;
286895fcc3eSMichael Lange           else numFaceVertices = s.nd;
287f7320561SMichael Lange         } else {              /* Data section */
288895fcc3eSMichael Lange           if (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices) {
289895fcc3eSMichael Lange             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported");
290f7320561SMichael Lange           }
291895fcc3eSMichael Lange           if (numFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file");
292f7320561SMichael Lange           if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd;
293f7320561SMichael Lange           numFaceEntries = numFaceVertices + 2;
294f7320561SMichael Lange           if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);}
295ec78a56aSMichael Lange           if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);}
296f7320561SMichael Lange           ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr);
297ec78a56aSMichael Lange           /* Record the zoneID for each face set */
298ec78a56aSMichael Lange           for (f = s.first -1; f < s.last; f++) faceZoneIDs[f] = s.zoneID;
299f7320561SMichael Lange           ierr = PetscFree(s.data);CHKERRQ(ierr);
300f7320561SMichael Lange         }
3012f0bd6dcSMichael Lange       }
3022f0bd6dcSMichael Lange     } while (s.index >= 0);
3032f0bd6dcSMichael Lange   }
304f7320561SMichael Lange   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
305f7320561SMichael Lange   if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension");
306f7320561SMichael Lange 
307f7320561SMichael Lange   /* Allocate cell-vertex mesh */
308f7320561SMichael Lange   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
309f7320561SMichael Lange   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
310f7320561SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
311f7320561SMichael Lange   ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr);
312f7320561SMichael Lange   if (!rank) {
313895fcc3eSMichael Lange     if (numCells < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file");
314895fcc3eSMichael Lange     /* If no cell type was given we assume simplices */
315895fcc3eSMichael Lange     if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1;
316f7320561SMichael Lange     for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);}
317f7320561SMichael Lange   }
318f7320561SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
319f7320561SMichael Lange 
320f7320561SMichael Lange   if (!rank) {
321f7320561SMichael Lange     /* Derive cell-vertex list from face-vertex and face-cell maps */
322f7320561SMichael Lange     ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr);
323f7320561SMichael Lange     for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1;
324f7320561SMichael Lange     for (f = 0; f < numFaces; f++) {
325f7320561SMichael Lange       PetscInt *cell;
326f7320561SMichael Lange       const PetscInt cl = faces[f*numFaceEntries + numFaceVertices];
327f7320561SMichael Lange       const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1];
328f7320561SMichael Lange       const PetscInt *face = &(faces[f*numFaceEntries]);
329f7320561SMichael Lange 
330f7320561SMichael Lange       if (cl > 0) {
331f7320561SMichael Lange         cell = &(cellVertices[(cl-1) * numCellVertices]);
332f7320561SMichael Lange         for (v = 0; v < numFaceVertices; v++) {
333f7320561SMichael Lange           PetscBool found = PETSC_FALSE;
334f7320561SMichael Lange           for (c = 0; c < numCellVertices; c++) {
335f7320561SMichael Lange             if (cell[c] < 0) break;
336f7320561SMichael Lange             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
337f7320561SMichael Lange           }
338f7320561SMichael Lange           if (!found) cell[c] = face[v]-1 + numCells;
339f7320561SMichael Lange         }
340f7320561SMichael Lange       }
341f7320561SMichael Lange       if (cr > 0) {
342f7320561SMichael Lange         cell = &(cellVertices[(cr-1) * numCellVertices]);
343f7320561SMichael Lange         for (v = 0; v < numFaceVertices; v++) {
344f7320561SMichael Lange           PetscBool found = PETSC_FALSE;
345f7320561SMichael Lange           for (c = 0; c < numCellVertices; c++) {
346f7320561SMichael Lange             if (cell[c] < 0) break;
347f7320561SMichael Lange             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
348f7320561SMichael Lange           }
349f7320561SMichael Lange           if (!found) cell[c] = face[v]-1 + numCells;
350f7320561SMichael Lange         }
351f7320561SMichael Lange       }
352f7320561SMichael Lange     }
353f7320561SMichael Lange     for (c = 0; c < numCells; c++) {
354f7320561SMichael Lange       ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr);
355f7320561SMichael Lange     }
356f7320561SMichael Lange   }
357f7320561SMichael Lange   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
358f7320561SMichael Lange   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
359f7320561SMichael Lange   if (interpolate) {
360f7320561SMichael Lange     DM idm = NULL;
361f7320561SMichael Lange 
362f7320561SMichael Lange     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
363f7320561SMichael Lange     ierr = DMDestroy(dm);CHKERRQ(ierr);
364f7320561SMichael Lange     *dm  = idm;
365f7320561SMichael Lange   }
366f7320561SMichael Lange 
367ec78a56aSMichael Lange   if (!rank) {
368631eb916SMichael Lange     PetscInt fi, joinSize, meetSize, *fverts, cells[2];
369631eb916SMichael Lange     const PetscInt *join, *meet;
370ec78a56aSMichael Lange     ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr);
371ec78a56aSMichael Lange     /* Mark facets by finding the full join of all adjacent vertices */
372ec78a56aSMichael Lange     for (f = 0; f < numFaces; f++) {
373631eb916SMichael Lange       const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1;
374631eb916SMichael Lange       const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1;
375631eb916SMichael Lange       if (cl > 0 && cr > 0) {
376631eb916SMichael Lange         /* If we know both adjoining cells we can use a single-level meet */
377631eb916SMichael Lange         cells[0] = cl; cells[1] = cr;
378631eb916SMichael Lange         ierr = DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet);CHKERRQ(ierr);
379631eb916SMichael Lange         if (meetSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f);
380631eb916SMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", meet[0], faceZoneIDs[f]);CHKERRQ(ierr);
381631eb916SMichael Lange         ierr = DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet);CHKERRQ(ierr);
382631eb916SMichael Lange       } else {
383ec78a56aSMichael Lange         for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1;
384ec78a56aSMichael Lange         ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr);
385ec78a56aSMichael Lange         if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f);
386ec78a56aSMichael Lange         ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr);
387ec78a56aSMichael Lange         ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr);
388ec78a56aSMichael Lange       }
389631eb916SMichael Lange     }
390ec78a56aSMichael Lange     ierr = PetscFree(fverts);CHKERRQ(ierr);
391ec78a56aSMichael Lange   }
392ec78a56aSMichael Lange 
3933f6dc66eSMichael Lange   /* Read coordinates */
3943f6dc66eSMichael Lange   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
3953f6dc66eSMichael Lange   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
3963f6dc66eSMichael Lange   ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr);
3973f6dc66eSMichael Lange   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
3983f6dc66eSMichael Lange   for (v = numCells; v < numCells+numVertices; ++v) {
3993f6dc66eSMichael Lange     ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr);
4003f6dc66eSMichael Lange     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr);
4013f6dc66eSMichael Lange   }
4023f6dc66eSMichael Lange   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
4033f6dc66eSMichael Lange   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
4043f6dc66eSMichael Lange   ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr);
4053f6dc66eSMichael Lange   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
4063f6dc66eSMichael Lange   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
4073f6dc66eSMichael Lange   ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr);
4083f6dc66eSMichael Lange   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
4093f6dc66eSMichael Lange   if (!rank) {
4103f6dc66eSMichael Lange     for (v = 0; v < numVertices; ++v) {
4113f6dc66eSMichael Lange       for (d = 0; d < dim; ++d) {
4123f6dc66eSMichael Lange         coords[v*dim+d] = coordsIn[v*dim+d];
4133f6dc66eSMichael Lange       }
4143f6dc66eSMichael Lange     }
4153f6dc66eSMichael Lange   }
4163f6dc66eSMichael Lange   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
4173f6dc66eSMichael Lange   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
4183f6dc66eSMichael Lange   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
419f7320561SMichael Lange   if (!rank) {
420f7320561SMichael Lange     ierr = PetscFree(cellVertices);CHKERRQ(ierr);
421f7320561SMichael Lange     ierr = PetscFree(faces);CHKERRQ(ierr);
422ec78a56aSMichael Lange     ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr);
4233f6dc66eSMichael Lange     ierr = PetscFree(coordsIn);CHKERRQ(ierr);
424f7320561SMichael Lange   }
4252f0bd6dcSMichael Lange   PetscFunctionReturn(0);
4262f0bd6dcSMichael Lange }
427