10039db0dSBarry Smith #define PETSC_DESIRE_FEATURE_TEST_MACROS /* for fileno() */ 22f0bd6dcSMichael Lange #define PETSCDM_DLL 3af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 42f0bd6dcSMichael Lange 52f0bd6dcSMichael Lange /*@C 62f0bd6dcSMichael Lange DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file 72f0bd6dcSMichael Lange 82f0bd6dcSMichael Lange + comm - The MPI communicator 92f0bd6dcSMichael Lange . filename - Name of the Fluent mesh file 102f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 112f0bd6dcSMichael Lange 122f0bd6dcSMichael Lange Output Parameter: 132f0bd6dcSMichael Lange . dm - The DM object representing the mesh 142f0bd6dcSMichael Lange 152f0bd6dcSMichael Lange Level: beginner 162f0bd6dcSMichael Lange 172f0bd6dcSMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate() 182f0bd6dcSMichael Lange @*/ 192f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 202f0bd6dcSMichael Lange { 212f0bd6dcSMichael Lange PetscViewer viewer; 222f0bd6dcSMichael Lange PetscErrorCode ierr; 232f0bd6dcSMichael Lange 242f0bd6dcSMichael Lange PetscFunctionBegin; 252f0bd6dcSMichael Lange /* Create file viewer and build plex */ 262f0bd6dcSMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 272f0bd6dcSMichael Lange ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 282f0bd6dcSMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 292f0bd6dcSMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 302f0bd6dcSMichael Lange ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr); 312f0bd6dcSMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 322f0bd6dcSMichael Lange PetscFunctionReturn(0); 332f0bd6dcSMichael Lange } 342f0bd6dcSMichael Lange 359dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 362f0bd6dcSMichael Lange { 374b375a39SMichael Lange PetscInt ret, i = 0; 382f0bd6dcSMichael Lange PetscErrorCode ierr; 392f0bd6dcSMichael Lange 402f0bd6dcSMichael Lange PetscFunctionBegin; 414b375a39SMichael Lange do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR);CHKERRQ(ierr);} 424b375a39SMichael Lange while (ret > 0 && buffer[i-1] != '\0' && buffer[i-1] != delim); 4336342426SSatish Balay if (!ret) buffer[i-1] = '\0'; else buffer[i] = '\0'; 442f0bd6dcSMichael Lange PetscFunctionReturn(0); 452f0bd6dcSMichael Lange } 462f0bd6dcSMichael Lange 479dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary) 4819d58f9dSMichael Lange { 4970c9a859SSatish Balay int fdes=0; 5019d58f9dSMichael Lange FILE *file; 5119d58f9dSMichael Lange PetscInt i; 5219d58f9dSMichael Lange PetscErrorCode ierr; 5319d58f9dSMichael Lange 5419d58f9dSMichael Lange PetscFunctionBegin; 5519d58f9dSMichael Lange if (binary) { 5619d58f9dSMichael Lange /* Extract raw file descriptor to read binary block */ 5719d58f9dSMichael Lange ierr = PetscViewerASCIIGetPointer(viewer, &file);CHKERRQ(ierr); 5819d58f9dSMichael Lange fflush(file); fdes = fileno(file); 5919d58f9dSMichael Lange } 6019d58f9dSMichael Lange 6119d58f9dSMichael Lange if (!binary && dtype == PETSC_INT) { 6219d58f9dSMichael Lange char cbuf[256]; 63cfb60857SMatthew G. Knepley unsigned int ibuf; 64cfb60857SMatthew G. Knepley int snum; 6519d58f9dSMichael Lange /* Parse hexadecimal ascii integers */ 6619d58f9dSMichael Lange for (i = 0; i < count; i++) { 67060da220SMatthew G. Knepley ierr = PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING);CHKERRQ(ierr); 6819d58f9dSMichael Lange snum = sscanf(cbuf, "%x", &ibuf); 69*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 7019d58f9dSMichael Lange ((PetscInt*)data)[i] = (PetscInt)ibuf; 7119d58f9dSMichael Lange } 7219d58f9dSMichael Lange } else if (binary && dtype == PETSC_INT) { 7319d58f9dSMichael Lange /* Always read 32-bit ints and cast to PetscInt */ 7419d58f9dSMichael Lange int *ibuf; 7519d58f9dSMichael Lange ierr = PetscMalloc1(count, &ibuf);CHKERRQ(ierr); 769860990eSLisandro Dalcin ierr = PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM);CHKERRQ(ierr); 7719d58f9dSMichael Lange ierr = PetscByteSwap(ibuf, PETSC_ENUM, count);CHKERRQ(ierr); 7819d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]); 7919d58f9dSMichael Lange ierr = PetscFree(ibuf);CHKERRQ(ierr); 8019d58f9dSMichael Lange 8119d58f9dSMichael Lange } else if (binary && dtype == PETSC_SCALAR) { 8219d58f9dSMichael Lange float *fbuf; 8319d58f9dSMichael Lange /* Always read 32-bit floats and cast to PetscScalar */ 8419d58f9dSMichael Lange ierr = PetscMalloc1(count, &fbuf);CHKERRQ(ierr); 859860990eSLisandro Dalcin ierr = PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT);CHKERRQ(ierr); 8619d58f9dSMichael Lange ierr = PetscByteSwap(fbuf, PETSC_FLOAT, count);CHKERRQ(ierr); 8719d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]); 8819d58f9dSMichael Lange ierr = PetscFree(fbuf);CHKERRQ(ierr); 8919d58f9dSMichael Lange } else { 90060da220SMatthew G. Knepley ierr = PetscViewerASCIIRead(viewer, data, count, NULL, dtype);CHKERRQ(ierr); 9119d58f9dSMichael Lange } 9219d58f9dSMichael Lange PetscFunctionReturn(0); 9319d58f9dSMichael Lange } 9419d58f9dSMichael Lange 959dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 962f0bd6dcSMichael Lange { 972f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 982f0bd6dcSMichael Lange int snum; 992f0bd6dcSMichael Lange PetscErrorCode ierr; 1002f0bd6dcSMichael Lange 1012f0bd6dcSMichael Lange PetscFunctionBegin; 1022f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 1032f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 1042f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr); 1052f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->index)); 1062f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 1072f0bd6dcSMichael Lange if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} 1082f0bd6dcSMichael Lange 1092f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 1102f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1112f0bd6dcSMichael Lange 1122f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 1132f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1142f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->nd)); 115*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1162f0bd6dcSMichael Lange 11719d58f9dSMichael Lange } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 1182f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1192f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 120*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1213f6dc66eSMichael Lange if (s->zoneID > 0) { 1223f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 1233f6dc66eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 1243f6dc66eSMichael Lange ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); 1255c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 12619d58f9dSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1273f6dc66eSMichael Lange } 12819d58f9dSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1292f0bd6dcSMichael Lange 13019d58f9dSMichael Lange } else if (s->index == 12 || s->index == 2012) { /* Cells */ 1312f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1322f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 133*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1342f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1352f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 136*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 4,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 137f7320561SMichael Lange } else { /* Data section */ 1382f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 139*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 140895fcc3eSMichael Lange if (s->nd == 0) { 141895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 14219d58f9dSMichael Lange PetscInt numCells = s->last - s->first + 1; 143895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 14419d58f9dSMichael Lange ierr = PetscMalloc1(numCells, (PetscInt**)&s->data);CHKERRQ(ierr); 1455c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 14619d58f9dSMichael Lange ierr = PetscFree(s->data);CHKERRQ(ierr); 147895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 148895fcc3eSMichael Lange } 1492f0bd6dcSMichael Lange } 1502f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1512f0bd6dcSMichael Lange 15219d58f9dSMichael Lange } else if (s->index == 13 || s->index == 2013) { /* Faces */ 1532f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1542f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 155*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1562f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1572f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 158*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 4,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 159f7320561SMichael Lange } else { /* Data section */ 160137d0469SJed Brown PetscInt f, numEntries, numFaces; 1612f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 162*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 163f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 164f7320561SMichael Lange switch (s->nd) { 165895fcc3eSMichael Lange case 0: numEntries = PETSC_DETERMINE; break; 166f7320561SMichael Lange case 2: numEntries = 2 + 2; break; /* linear */ 167f7320561SMichael Lange case 3: numEntries = 2 + 3; break; /* triangular */ 168f7320561SMichael Lange case 4: numEntries = 2 + 4; break; /* quadrilateral */ 169f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 170f7320561SMichael Lange } 171f7320561SMichael Lange numFaces = s->last-s->first + 1; 172895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 173895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 174f7320561SMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 175895fcc3eSMichael Lange } 176f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 177895fcc3eSMichael Lange if (s->nd == 0) { 178895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 179137d0469SJed Brown PetscInt numFaceVert = 0; 1805c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 181895fcc3eSMichael Lange if (numEntries == PETSC_DETERMINE) { 182895fcc3eSMichael Lange numEntries = numFaceVert + 2; 183895fcc3eSMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 184895fcc3eSMichael Lange } else { 185*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numEntries != numFaceVert + 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 186895fcc3eSMichael Lange } 187895fcc3eSMichael Lange } 1885c0a1baaSMatthew G. Knepley ierr = DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE);CHKERRQ(ierr); 189f7320561SMichael Lange } 190895fcc3eSMichael Lange s->nd = numEntries - 2; 191f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1922f0bd6dcSMichael Lange } 1932f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1942f0bd6dcSMichael Lange 1952f0bd6dcSMichael Lange } else { /* Unknown section type */ 196060da220SMatthew G. Knepley PetscInt depth = 1; 1972f0bd6dcSMichael Lange do { 1982f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 199060da220SMatthew G. Knepley do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR);CHKERRQ(ierr);} 2002f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 2012f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 2022f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 2032f0bd6dcSMichael Lange } while (depth > 0); 2042f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); 2052f0bd6dcSMichael Lange } 2062f0bd6dcSMichael Lange PetscFunctionReturn(0); 2072f0bd6dcSMichael Lange } 2082f0bd6dcSMichael Lange 2092f0bd6dcSMichael Lange /*@C 2102f0bd6dcSMichael Lange DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 2112f0bd6dcSMichael Lange 212d083f849SBarry Smith Collective 2132f0bd6dcSMichael Lange 2142f0bd6dcSMichael Lange Input Parameters: 2152f0bd6dcSMichael Lange + comm - The MPI communicator 2162f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file 2172f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 2182f0bd6dcSMichael Lange 2192f0bd6dcSMichael Lange Output Parameter: 2202f0bd6dcSMichael Lange . dm - The DM object representing the mesh 2212f0bd6dcSMichael Lange 2222f0bd6dcSMichael Lange Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 2232f0bd6dcSMichael Lange 2242f0bd6dcSMichael Lange Level: beginner 2252f0bd6dcSMichael Lange 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; 2347368db69SLisandro Dalcin DMLabel faceSets = NULL; 2353f6dc66eSMichael Lange PetscInt d, coordSize; 2363f6dc66eSMichael Lange PetscScalar *coords, *coordsIn = NULL; 2373f6dc66eSMichael Lange PetscSection coordSection; 2383f6dc66eSMichael Lange Vec coordinates; 2392f0bd6dcSMichael Lange PetscErrorCode ierr; 2402f0bd6dcSMichael Lange 2412f0bd6dcSMichael Lange PetscFunctionBegin; 242ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 2432f0bd6dcSMichael Lange 244dd400576SPatrick Sanan if (rank == 0) { 2452f0bd6dcSMichael Lange FluentSection s; 2462f0bd6dcSMichael Lange do { 2472f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr); 2482f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 249f7320561SMichael Lange dim = s.nd; 2502f0bd6dcSMichael Lange 25119d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 252f7320561SMichael Lange if (s.zoneID == 0) numVertices = s.last; 2533f6dc66eSMichael Lange else { 254*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(coordsIn,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 2551a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *) s.data; 2563f6dc66eSMichael Lange } 2572f0bd6dcSMichael Lange 25819d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 259f7320561SMichael Lange if (s.zoneID == 0) numCells = s.last; 260f7320561SMichael Lange else { 261f7320561SMichael Lange switch (s.nd) { 262895fcc3eSMichael Lange case 0: numCellVertices = PETSC_DETERMINE; break; 263f7320561SMichael Lange case 1: numCellVertices = 3; break; /* triangular */ 264f7320561SMichael Lange case 2: numCellVertices = 4; break; /* tetrahedral */ 265f7320561SMichael Lange case 3: numCellVertices = 4; break; /* quadrilateral */ 266f7320561SMichael Lange case 4: numCellVertices = 8; break; /* hexahedral */ 267f7320561SMichael Lange case 5: numCellVertices = 5; break; /* pyramid */ 268f7320561SMichael Lange case 6: numCellVertices = 6; break; /* wedge */ 269895fcc3eSMichael Lange default: numCellVertices = PETSC_DETERMINE; 270f7320561SMichael Lange } 271f7320561SMichael Lange } 2722f0bd6dcSMichael Lange 27319d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 274f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 275930bae4bSMatthew G. Knepley numFaces = (PetscInt) (s.last - s.first + 1); 27635462f7fSMichael Lange if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 277895fcc3eSMichael Lange else numFaceVertices = s.nd; 278f7320561SMichael Lange } else { /* Data section */ 279cf554e2cSMatthew G. Knepley unsigned int z; 280cf554e2cSMatthew G. Knepley 281*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 282*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numFaces < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 283f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 284f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 285f7320561SMichael Lange if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);} 286ec78a56aSMichael Lange if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);} 287580bdb30SBarry Smith ierr = PetscMemcpy(&faces[(s.first-1)*numFaceEntries], s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr); 288ec78a56aSMichael Lange /* Record the zoneID for each face set */ 289cf554e2cSMatthew G. Knepley for (z = s.first -1; z < s.last; z++) faceZoneIDs[z] = s.zoneID; 290f7320561SMichael Lange ierr = PetscFree(s.data);CHKERRQ(ierr); 291f7320561SMichael Lange } 2922f0bd6dcSMichael Lange } 2932f0bd6dcSMichael Lange } while (s.index >= 0); 2942f0bd6dcSMichael Lange } 295ffc4695bSBarry Smith ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRMPI(ierr); 296*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 297f7320561SMichael Lange 298f7320561SMichael Lange /* Allocate cell-vertex mesh */ 299f7320561SMichael Lange ierr = DMCreate(comm, dm);CHKERRQ(ierr); 300f7320561SMichael Lange ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 301f7320561SMichael Lange ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 302f7320561SMichael Lange ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr); 303dd400576SPatrick Sanan if (rank == 0) { 304*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numCells < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 305895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 306895fcc3eSMichael Lange if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1; 307f7320561SMichael Lange for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);} 308f7320561SMichael Lange } 309f7320561SMichael Lange ierr = DMSetUp(*dm);CHKERRQ(ierr); 310f7320561SMichael Lange 311dd400576SPatrick Sanan if (rank == 0 && faces) { 312f7320561SMichael Lange /* Derive cell-vertex list from face-vertex and face-cell maps */ 313f7320561SMichael Lange ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr); 314f7320561SMichael Lange for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 315f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 316f7320561SMichael Lange PetscInt *cell; 317f7320561SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 318f7320561SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 319f7320561SMichael Lange const PetscInt *face = &(faces[f*numFaceEntries]); 320f7320561SMichael Lange 321f7320561SMichael Lange if (cl > 0) { 322f7320561SMichael Lange cell = &(cellVertices[(cl-1) * numCellVertices]); 323f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 324f7320561SMichael Lange PetscBool found = PETSC_FALSE; 325f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 326f7320561SMichael Lange if (cell[c] < 0) break; 327f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 328f7320561SMichael Lange } 329f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 330f7320561SMichael Lange } 331f7320561SMichael Lange } 332f7320561SMichael Lange if (cr > 0) { 333f7320561SMichael Lange cell = &(cellVertices[(cr-1) * numCellVertices]); 334f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 335f7320561SMichael Lange PetscBool found = PETSC_FALSE; 336f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 337f7320561SMichael Lange if (cell[c] < 0) break; 338f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 339f7320561SMichael Lange } 340f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 341f7320561SMichael Lange } 342f7320561SMichael Lange } 343f7320561SMichael Lange } 344f7320561SMichael Lange for (c = 0; c < numCells; c++) { 345f7320561SMichael Lange ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr); 346f7320561SMichael Lange } 347f7320561SMichael Lange } 348f7320561SMichael Lange ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 349f7320561SMichael Lange ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 350f7320561SMichael Lange if (interpolate) { 3515fd9971aSMatthew G. Knepley DM idm; 352f7320561SMichael Lange 353f7320561SMichael Lange ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 354f7320561SMichael Lange ierr = DMDestroy(dm);CHKERRQ(ierr); 355f7320561SMichael Lange *dm = idm; 356f7320561SMichael Lange } 357f7320561SMichael Lange 358dd400576SPatrick Sanan if (rank == 0 && faces) { 359631eb916SMichael Lange PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 360631eb916SMichael Lange const PetscInt *join, *meet; 361ec78a56aSMichael Lange ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr); 362ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 363ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 364631eb916SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1; 365631eb916SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1; 366631eb916SMichael Lange if (cl > 0 && cr > 0) { 367631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 368631eb916SMichael Lange cells[0] = cl; cells[1] = cr; 369631eb916SMichael Lange ierr = DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet);CHKERRQ(ierr); 370*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(meetSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 3717368db69SLisandro Dalcin ierr = DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], faceZoneIDs[f]);CHKERRQ(ierr); 372631eb916SMichael Lange ierr = DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet);CHKERRQ(ierr); 373631eb916SMichael Lange } else { 374ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1; 375ec78a56aSMichael Lange ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 376*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(joinSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 3777368db69SLisandro Dalcin ierr = DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr); 378ec78a56aSMichael Lange ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 379ec78a56aSMichael Lange } 380631eb916SMichael Lange } 381ec78a56aSMichael Lange ierr = PetscFree(fverts);CHKERRQ(ierr); 382ec78a56aSMichael Lange } 383ec78a56aSMichael Lange 3847368db69SLisandro Dalcin { /* Create Face Sets label at all processes */ 3857368db69SLisandro Dalcin enum {n = 1}; 3867368db69SLisandro Dalcin PetscBool flag[n]; 3877368db69SLisandro Dalcin 3887368db69SLisandro Dalcin flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE; 3897368db69SLisandro Dalcin ierr = MPI_Bcast(flag, n, MPIU_BOOL, 0, comm);CHKERRMPI(ierr); 3907368db69SLisandro Dalcin if (flag[0]) {ierr = DMCreateLabel(*dm, "Face Sets");CHKERRQ(ierr);} 3917368db69SLisandro Dalcin } 3927368db69SLisandro Dalcin 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); 4048b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &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); 409dd400576SPatrick Sanan if (rank == 0 && coordsIn) { 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); 4197368db69SLisandro Dalcin 420dd400576SPatrick Sanan if (rank == 0) { 4217368db69SLisandro Dalcin ierr = PetscFree(cellVertices);CHKERRQ(ierr); 422f7320561SMichael Lange ierr = PetscFree(faces);CHKERRQ(ierr); 423ec78a56aSMichael Lange ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr); 4243f6dc66eSMichael Lange ierr = PetscFree(coordsIn);CHKERRQ(ierr); 425f7320561SMichael Lange } 4262f0bd6dcSMichael Lange PetscFunctionReturn(0); 4272f0bd6dcSMichael Lange } 428