12f0bd6dcSMichael Lange #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 32f0bd6dcSMichael Lange 42f0bd6dcSMichael Lange /*@C 52f0bd6dcSMichael Lange DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file 62f0bd6dcSMichael Lange 72f0bd6dcSMichael Lange + comm - The MPI communicator 82f0bd6dcSMichael Lange . filename - Name of the Fluent mesh file 92f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 102f0bd6dcSMichael Lange 112f0bd6dcSMichael Lange Output Parameter: 122f0bd6dcSMichael Lange . dm - The DM object representing the mesh 132f0bd6dcSMichael Lange 142f0bd6dcSMichael Lange Level: beginner 152f0bd6dcSMichael Lange 162f0bd6dcSMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate() 172f0bd6dcSMichael Lange @*/ 182f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 192f0bd6dcSMichael Lange { 202f0bd6dcSMichael Lange PetscViewer viewer; 212f0bd6dcSMichael Lange PetscErrorCode ierr; 222f0bd6dcSMichael Lange 232f0bd6dcSMichael Lange PetscFunctionBegin; 242f0bd6dcSMichael Lange /* Create file viewer and build plex */ 252f0bd6dcSMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 262f0bd6dcSMichael Lange ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 272f0bd6dcSMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 282f0bd6dcSMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 292f0bd6dcSMichael Lange ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr); 302f0bd6dcSMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 312f0bd6dcSMichael Lange PetscFunctionReturn(0); 322f0bd6dcSMichael Lange } 332f0bd6dcSMichael Lange 34*9dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 352f0bd6dcSMichael Lange { 364b375a39SMichael Lange PetscInt ret, i = 0; 372f0bd6dcSMichael Lange PetscErrorCode ierr; 382f0bd6dcSMichael Lange 392f0bd6dcSMichael Lange PetscFunctionBegin; 404b375a39SMichael Lange do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR);CHKERRQ(ierr);} 414b375a39SMichael Lange while (ret > 0 && buffer[i-1] != '\0' && buffer[i-1] != delim); 422f0bd6dcSMichael Lange buffer[i] = '\0'; 432f0bd6dcSMichael Lange PetscFunctionReturn(0); 442f0bd6dcSMichael Lange } 452f0bd6dcSMichael Lange 46*9dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary) 4719d58f9dSMichael Lange { 4870c9a859SSatish Balay int fdes=0; 4919d58f9dSMichael Lange FILE *file; 5019d58f9dSMichael Lange PetscInt i; 5119d58f9dSMichael Lange PetscErrorCode ierr; 5219d58f9dSMichael Lange 5319d58f9dSMichael Lange PetscFunctionBegin; 5419d58f9dSMichael Lange if (binary) { 5519d58f9dSMichael Lange /* Extract raw file descriptor to read binary block */ 5619d58f9dSMichael Lange ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr); 5719d58f9dSMichael Lange fflush(file); fdes = fileno(file); 5819d58f9dSMichael Lange } 5919d58f9dSMichael Lange 6019d58f9dSMichael Lange if (!binary && dtype == PETSC_INT) { 6119d58f9dSMichael Lange char cbuf[256]; 62cfb60857SMatthew G. Knepley unsigned int ibuf; 63cfb60857SMatthew G. Knepley int snum; 6419d58f9dSMichael Lange /* Parse hexadecimal ascii integers */ 6519d58f9dSMichael Lange for (i = 0; i < count; i++) { 66060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 6719d58f9dSMichael Lange snum = sscanf(cbuf, "%x", &ibuf); 6819d58f9dSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 6919d58f9dSMichael Lange ((PetscInt*)data)[i] = (PetscInt)ibuf; 7019d58f9dSMichael Lange } 7119d58f9dSMichael Lange } else if (binary && dtype == PETSC_INT) { 7219d58f9dSMichael Lange /* Always read 32-bit ints and cast to PetscInt */ 7319d58f9dSMichael Lange int *ibuf; 7419d58f9dSMichael Lange ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr); 7519d58f9dSMichael Lange ierr = PetscBinaryRead(fdes, ibuf, count, PETSC_ENUM);CHKERRQ(ierr); 7619d58f9dSMichael Lange ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr); 7719d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]); 7819d58f9dSMichael Lange ierr = PetscFree(ibuf);CHKERRQ(ierr); 7919d58f9dSMichael Lange 8019d58f9dSMichael Lange } else if (binary && dtype == PETSC_SCALAR) { 8119d58f9dSMichael Lange float *fbuf; 8219d58f9dSMichael Lange /* Always read 32-bit floats and cast to PetscScalar */ 8319d58f9dSMichael Lange ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr); 8419d58f9dSMichael Lange ierr = PetscBinaryRead(fdes, fbuf, count, PETSC_FLOAT);CHKERRQ(ierr); 8519d58f9dSMichael Lange ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr); 8619d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]); 8719d58f9dSMichael Lange ierr = PetscFree(fbuf);CHKERRQ(ierr); 8819d58f9dSMichael Lange } else { 89060da220SMatthew G. Knepley ierr = PetscViewerASCIIRead(viewer, data, count, NULL, dtype);CHKERRQ(ierr); 9019d58f9dSMichael Lange } 9119d58f9dSMichael Lange PetscFunctionReturn(0); 9219d58f9dSMichael Lange } 9319d58f9dSMichael Lange 94*9dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 952f0bd6dcSMichael Lange { 962f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 972f0bd6dcSMichael Lange int snum; 982f0bd6dcSMichael Lange PetscErrorCode ierr; 992f0bd6dcSMichael Lange 1002f0bd6dcSMichael Lange PetscFunctionBegin; 1012f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 1022f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 1032f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr); 1042f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->index)); 1052f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 1062f0bd6dcSMichael Lange if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} 1072f0bd6dcSMichael Lange 1082f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 1092f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1102f0bd6dcSMichael Lange 1112f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 1122f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1132f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->nd)); 1142f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1152f0bd6dcSMichael Lange 11619d58f9dSMichael Lange } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 1172f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1182f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1192f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1203f6dc66eSMichael Lange if (s->zoneID > 0) { 1213f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 1223f6dc66eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 1233f6dc66eSMichael Lange ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); 1245c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 12519d58f9dSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1263f6dc66eSMichael Lange } 12719d58f9dSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1282f0bd6dcSMichael Lange 12919d58f9dSMichael Lange } else if (s->index == 12 || s->index == 2012) { /* Cells */ 1302f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1312f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1322f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1332f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1342f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1352f0bd6dcSMichael Lange if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 136f7320561SMichael Lange } else { /* Data section */ 1372f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1382f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 139895fcc3eSMichael Lange if (s->nd == 0) { 140895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 14119d58f9dSMichael Lange PetscInt numCells = s->last - s->first + 1; 142895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 14319d58f9dSMichael Lange ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr); 1445c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 14519d58f9dSMichael Lange ierr = PetscFree(s->data);CHKERRQ(ierr); 146895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 147895fcc3eSMichael Lange } 1482f0bd6dcSMichael Lange } 1492f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1502f0bd6dcSMichael Lange 15119d58f9dSMichael Lange } else if (s->index == 13 || s->index == 2013) { /* Faces */ 1522f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1532f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1542f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1552f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1562f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1572f0bd6dcSMichael Lange if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 158f7320561SMichael Lange } else { /* Data section */ 159137d0469SJed Brown PetscInt f, numEntries, numFaces; 1602f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1612f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 162f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 163f7320561SMichael Lange switch (s->nd) { 164895fcc3eSMichael Lange case 0: numEntries = PETSC_DETERMINE; break; 165f7320561SMichael Lange case 2: numEntries = 2 + 2; break; /* linear */ 166f7320561SMichael Lange case 3: numEntries = 2 + 3; break; /* triangular */ 167f7320561SMichael Lange case 4: numEntries = 2 + 4; break; /* quadrilateral */ 168f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 169f7320561SMichael Lange } 170f7320561SMichael Lange numFaces = s->last-s->first + 1; 171895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 172895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 173f7320561SMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 174895fcc3eSMichael Lange } 175f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 176895fcc3eSMichael Lange if (s->nd == 0) { 177895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 178137d0469SJed Brown PetscInt numFaceVert = 0; 1795c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 180895fcc3eSMichael Lange if (numEntries == PETSC_DETERMINE) { 181895fcc3eSMichael Lange numEntries = numFaceVert + 2; 182895fcc3eSMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 183895fcc3eSMichael Lange } else { 1846c4ed002SBarry Smith if (numEntries != numFaceVert + 2) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 185895fcc3eSMichael Lange } 186895fcc3eSMichael Lange } 1875c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 188f7320561SMichael Lange } 189895fcc3eSMichael Lange s->nd = numEntries - 2; 190f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1912f0bd6dcSMichael Lange } 1922f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1932f0bd6dcSMichael Lange 1942f0bd6dcSMichael Lange } else { /* Unknown section type */ 195060da220SMatthew G. Knepley PetscInt depth = 1; 1962f0bd6dcSMichael Lange do { 1972f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 198060da220SMatthew G. Knepley do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);} 1992f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 2002f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 2012f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 2022f0bd6dcSMichael Lange } while (depth > 0); 2032f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); 2042f0bd6dcSMichael Lange } 2052f0bd6dcSMichael Lange PetscFunctionReturn(0); 2062f0bd6dcSMichael Lange } 2072f0bd6dcSMichael Lange 2082f0bd6dcSMichael Lange /*@C 2092f0bd6dcSMichael Lange DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 2102f0bd6dcSMichael Lange 2112f0bd6dcSMichael Lange Collective on comm 2122f0bd6dcSMichael Lange 2132f0bd6dcSMichael Lange Input Parameters: 2142f0bd6dcSMichael Lange + comm - The MPI communicator 2152f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file 2162f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 2172f0bd6dcSMichael Lange 2182f0bd6dcSMichael Lange Output Parameter: 2192f0bd6dcSMichael Lange . dm - The DM object representing the mesh 2202f0bd6dcSMichael Lange 2212f0bd6dcSMichael Lange Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 2222f0bd6dcSMichael Lange 2232f0bd6dcSMichael Lange Level: beginner 2242f0bd6dcSMichael Lange 2252f0bd6dcSMichael Lange .keywords: mesh, fluent, case 2262f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate() 2272f0bd6dcSMichael Lange @*/ 2282f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 2292f0bd6dcSMichael Lange { 2302f0bd6dcSMichael Lange PetscMPIInt rank; 231c2fc583bSSatish Balay PetscInt c, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE; 232cf554e2cSMatthew G. Knepley PetscInt numFaces = PETSC_DETERMINE, f, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE; 2333ed8799eSMichael Lange PetscInt *faces = NULL, *cellVertices = NULL, *faceZoneIDs = NULL; 2343f6dc66eSMichael Lange PetscInt d, coordSize; 2353f6dc66eSMichael Lange PetscScalar *coords, *coordsIn = NULL; 2363f6dc66eSMichael Lange PetscSection coordSection; 2373f6dc66eSMichael Lange Vec coordinates; 2382f0bd6dcSMichael Lange PetscErrorCode ierr; 2392f0bd6dcSMichael Lange 2402f0bd6dcSMichael Lange PetscFunctionBegin; 2412f0bd6dcSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2422f0bd6dcSMichael Lange 2432f0bd6dcSMichael Lange if (!rank) { 2442f0bd6dcSMichael Lange FluentSection s; 2452f0bd6dcSMichael Lange do { 2462f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr); 2472f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 248f7320561SMichael Lange dim = s.nd; 2492f0bd6dcSMichael Lange 25019d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 251f7320561SMichael Lange if (s.zoneID == 0) numVertices = s.last; 2523f6dc66eSMichael Lange else { 2533f6dc66eSMichael Lange if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 2541a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *) s.data; 2553f6dc66eSMichael Lange } 2562f0bd6dcSMichael Lange 25719d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 258f7320561SMichael Lange if (s.zoneID == 0) numCells = s.last; 259f7320561SMichael Lange else { 260f7320561SMichael Lange switch (s.nd) { 261895fcc3eSMichael Lange case 0: numCellVertices = PETSC_DETERMINE; break; 262f7320561SMichael Lange case 1: numCellVertices = 3; break; /* triangular */ 263f7320561SMichael Lange case 2: numCellVertices = 4; break; /* tetrahedral */ 264f7320561SMichael Lange case 3: numCellVertices = 4; break; /* quadrilateral */ 265f7320561SMichael Lange case 4: numCellVertices = 8; break; /* hexahedral */ 266f7320561SMichael Lange case 5: numCellVertices = 5; break; /* pyramid */ 267f7320561SMichael Lange case 6: numCellVertices = 6; break; /* wedge */ 268895fcc3eSMichael Lange default: numCellVertices = PETSC_DETERMINE; 269f7320561SMichael Lange } 270f7320561SMichael Lange } 2712f0bd6dcSMichael Lange 27219d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 273f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 274930bae4bSMatthew G. Knepley numFaces = (PetscInt) (s.last - s.first + 1); 27535462f7fSMichael Lange if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 276895fcc3eSMichael Lange else numFaceVertices = s.nd; 277f7320561SMichael Lange } else { /* Data section */ 278cf554e2cSMatthew G. Knepley unsigned int z; 279cf554e2cSMatthew G. Knepley 2806c4ed002SBarry Smith if (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 281895fcc3eSMichael Lange if (numFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 282f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 283f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 284f7320561SMichael Lange if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);} 285ec78a56aSMichael Lange if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);} 286f7320561SMichael Lange ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr); 287ec78a56aSMichael Lange /* Record the zoneID for each face set */ 288cf554e2cSMatthew G. Knepley for (z = s.first -1; z < s.last; z++) faceZoneIDs[z] = s.zoneID; 289f7320561SMichael Lange ierr = PetscFree(s.data);CHKERRQ(ierr); 290f7320561SMichael Lange } 2912f0bd6dcSMichael Lange } 2922f0bd6dcSMichael Lange } while (s.index >= 0); 2932f0bd6dcSMichael Lange } 294f7320561SMichael Lange ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 295f7320561SMichael Lange if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 296f7320561SMichael Lange 297f7320561SMichael Lange /* Allocate cell-vertex mesh */ 298f7320561SMichael Lange ierr = DMCreate(comm, dm);CHKERRQ(ierr); 299f7320561SMichael Lange ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 300f7320561SMichael Lange ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 301f7320561SMichael Lange ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr); 302f7320561SMichael Lange if (!rank) { 303895fcc3eSMichael Lange if (numCells < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 304895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 305895fcc3eSMichael Lange if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1; 306f7320561SMichael Lange for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);} 307f7320561SMichael Lange } 308f7320561SMichael Lange ierr = DMSetUp(*dm);CHKERRQ(ierr); 309f7320561SMichael Lange 3103ed8799eSMichael Lange if (!rank && faces) { 311f7320561SMichael Lange /* Derive cell-vertex list from face-vertex and face-cell maps */ 312f7320561SMichael Lange ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr); 313f7320561SMichael Lange for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 314f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 315f7320561SMichael Lange PetscInt *cell; 316f7320561SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 317f7320561SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 318f7320561SMichael Lange const PetscInt *face = &(faces[f*numFaceEntries]); 319f7320561SMichael Lange 320f7320561SMichael Lange if (cl > 0) { 321f7320561SMichael Lange cell = &(cellVertices[(cl-1) * numCellVertices]); 322f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 323f7320561SMichael Lange PetscBool found = PETSC_FALSE; 324f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 325f7320561SMichael Lange if (cell[c] < 0) break; 326f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 327f7320561SMichael Lange } 328f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 329f7320561SMichael Lange } 330f7320561SMichael Lange } 331f7320561SMichael Lange if (cr > 0) { 332f7320561SMichael Lange cell = &(cellVertices[(cr-1) * numCellVertices]); 333f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 334f7320561SMichael Lange PetscBool found = PETSC_FALSE; 335f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 336f7320561SMichael Lange if (cell[c] < 0) break; 337f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 338f7320561SMichael Lange } 339f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 340f7320561SMichael Lange } 341f7320561SMichael Lange } 342f7320561SMichael Lange } 343f7320561SMichael Lange for (c = 0; c < numCells; c++) { 344f7320561SMichael Lange ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr); 345f7320561SMichael Lange } 346f7320561SMichael Lange } 347f7320561SMichael Lange ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 348f7320561SMichael Lange ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 349f7320561SMichael Lange if (interpolate) { 350f7320561SMichael Lange DM idm = NULL; 351f7320561SMichael Lange 352f7320561SMichael Lange ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 353f7320561SMichael Lange ierr = DMDestroy(dm);CHKERRQ(ierr); 354f7320561SMichael Lange *dm = idm; 355f7320561SMichael Lange } 356f7320561SMichael Lange 3573ed8799eSMichael Lange if (!rank && faces) { 358631eb916SMichael Lange PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 359631eb916SMichael Lange const PetscInt *join, *meet; 360ec78a56aSMichael Lange ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr); 361ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 362ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 363631eb916SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1; 364631eb916SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1; 365631eb916SMichael Lange if (cl > 0 && cr > 0) { 366631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 367631eb916SMichael Lange cells[0] = cl; cells[1] = cr; 368631eb916SMichael Lange ierr = DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet);CHKERRQ(ierr); 369631eb916SMichael Lange if (meetSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 370c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", meet[0], faceZoneIDs[f]);CHKERRQ(ierr); 371631eb916SMichael Lange ierr = DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet);CHKERRQ(ierr); 372631eb916SMichael Lange } else { 373ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1; 374ec78a56aSMichael Lange ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 375ec78a56aSMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 376c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr); 377ec78a56aSMichael Lange ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 378ec78a56aSMichael Lange } 379631eb916SMichael Lange } 380ec78a56aSMichael Lange ierr = PetscFree(fverts);CHKERRQ(ierr); 381ec78a56aSMichael Lange } 382ec78a56aSMichael Lange 3833f6dc66eSMichael Lange /* Read coordinates */ 3843f6dc66eSMichael Lange ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 3853f6dc66eSMichael Lange ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 3863f6dc66eSMichael Lange ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 3873f6dc66eSMichael Lange ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 3883f6dc66eSMichael Lange for (v = numCells; v < numCells+numVertices; ++v) { 3893f6dc66eSMichael Lange ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 3903f6dc66eSMichael Lange ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 3913f6dc66eSMichael Lange } 3923f6dc66eSMichael Lange ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 3933f6dc66eSMichael Lange ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 3948b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 3953f6dc66eSMichael Lange ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 3963f6dc66eSMichael Lange ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 3973f6dc66eSMichael Lange ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 3983f6dc66eSMichael Lange ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 3993ed8799eSMichael Lange if (!rank && coordsIn) { 4003f6dc66eSMichael Lange for (v = 0; v < numVertices; ++v) { 4013f6dc66eSMichael Lange for (d = 0; d < dim; ++d) { 4023f6dc66eSMichael Lange coords[v*dim+d] = coordsIn[v*dim+d]; 4033f6dc66eSMichael Lange } 4043f6dc66eSMichael Lange } 4053f6dc66eSMichael Lange } 4063f6dc66eSMichael Lange ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 4073f6dc66eSMichael Lange ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 4083f6dc66eSMichael Lange ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 409f7320561SMichael Lange if (!rank) { 4103ed8799eSMichael Lange if (cellVertices) {ierr = PetscFree(cellVertices);CHKERRQ(ierr);} 411f7320561SMichael Lange ierr = PetscFree(faces);CHKERRQ(ierr); 412ec78a56aSMichael Lange ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr); 4133f6dc66eSMichael Lange ierr = PetscFree(coordsIn);CHKERRQ(ierr); 414f7320561SMichael Lange } 4152f0bd6dcSMichael Lange PetscFunctionReturn(0); 4162f0bd6dcSMichael Lange } 417