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 232f0bd6dcSMichael Lange PetscFunctionBegin; 242f0bd6dcSMichael Lange /* Create file viewer and build plex */ 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerCreate(comm, &viewer)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFileSetName(viewer, filename)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent(comm, viewer, interpolate, dm)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 312f0bd6dcSMichael Lange PetscFunctionReturn(0); 322f0bd6dcSMichael Lange } 332f0bd6dcSMichael Lange 349dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 352f0bd6dcSMichael Lange { 364b375a39SMichael Lange PetscInt ret, i = 0; 372f0bd6dcSMichael Lange 382f0bd6dcSMichael Lange PetscFunctionBegin; 395f80ce2aSJacob Faibussowitsch do CHKERRQ(PetscViewerRead(viewer, &(buffer[i++]), 1, &ret, PETSC_CHAR)); 404b375a39SMichael Lange while (ret > 0 && buffer[i-1] != '\0' && buffer[i-1] != delim); 4136342426SSatish Balay if (!ret) buffer[i-1] = '\0'; else buffer[i] = '\0'; 422f0bd6dcSMichael Lange PetscFunctionReturn(0); 432f0bd6dcSMichael Lange } 442f0bd6dcSMichael Lange 459dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary) 4619d58f9dSMichael Lange { 4770c9a859SSatish Balay int fdes=0; 4819d58f9dSMichael Lange FILE *file; 4919d58f9dSMichael Lange PetscInt i; 5019d58f9dSMichael Lange 5119d58f9dSMichael Lange PetscFunctionBegin; 5219d58f9dSMichael Lange if (binary) { 5319d58f9dSMichael Lange /* Extract raw file descriptor to read binary block */ 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIGetPointer(viewer, &file)); 5519d58f9dSMichael Lange fflush(file); fdes = fileno(file); 5619d58f9dSMichael Lange } 5719d58f9dSMichael Lange 5819d58f9dSMichael Lange if (!binary && dtype == PETSC_INT) { 5919d58f9dSMichael Lange char cbuf[256]; 60cfb60857SMatthew G. Knepley unsigned int ibuf; 61cfb60857SMatthew G. Knepley int snum; 6219d58f9dSMichael Lange /* Parse hexadecimal ascii integers */ 6319d58f9dSMichael Lange for (i = 0; i < count; i++) { 645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING)); 6519d58f9dSMichael Lange snum = sscanf(cbuf, "%x", &ibuf); 662c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 6719d58f9dSMichael Lange ((PetscInt*)data)[i] = (PetscInt)ibuf; 6819d58f9dSMichael Lange } 6919d58f9dSMichael Lange } else if (binary && dtype == PETSC_INT) { 7019d58f9dSMichael Lange /* Always read 32-bit ints and cast to PetscInt */ 7119d58f9dSMichael Lange int *ibuf; 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(count, &ibuf)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscByteSwap(ibuf, PETSC_ENUM, count)); 7519d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscInt*)data)[i] = (PetscInt)(ibuf[i]); 765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(ibuf)); 7719d58f9dSMichael Lange 7819d58f9dSMichael Lange } else if (binary && dtype == PETSC_SCALAR) { 7919d58f9dSMichael Lange float *fbuf; 8019d58f9dSMichael Lange /* Always read 32-bit floats and cast to PetscScalar */ 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(count, &fbuf)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscByteSwap(fbuf, PETSC_FLOAT, count)); 8419d58f9dSMichael Lange for (i = 0; i < count; i++) ((PetscScalar*)data)[i] = (PetscScalar)(fbuf[i]); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(fbuf)); 8619d58f9dSMichael Lange } else { 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIRead(viewer, data, count, NULL, dtype)); 8819d58f9dSMichael Lange } 8919d58f9dSMichael Lange PetscFunctionReturn(0); 9019d58f9dSMichael Lange } 9119d58f9dSMichael Lange 929dad6f3fSMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 932f0bd6dcSMichael Lange { 942f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 952f0bd6dcSMichael Lange int snum; 962f0bd6dcSMichael Lange 972f0bd6dcSMichael Lange PetscFunctionBegin; 982f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 995f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ' ')); 1012f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->index)); 1022f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 1032f0bd6dcSMichael Lange if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} 1042f0bd6dcSMichael Lange 1052f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1072f0bd6dcSMichael Lange 1082f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1102f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->nd)); 1112c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1122f0bd6dcSMichael Lange 11319d58f9dSMichael Lange } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 1145f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1152f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1162c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1173f6dc66eSMichael Lange if (s->zoneID > 0) { 1183f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 1195f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd*numCoords, PETSC_SCALAR, s->index==2010 ? PETSC_TRUE : PETSC_FALSE)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1233f6dc66eSMichael Lange } 1245f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1252f0bd6dcSMichael Lange 12619d58f9dSMichael Lange } else if (s->index == 12 || s->index == 2012) { /* Cells */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1282f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1292c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1302f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1312f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1322c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 4,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 133f7320561SMichael Lange } else { /* Data section */ 1342f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1352c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 136895fcc3eSMichael Lange if (s->nd == 0) { 137895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 13819d58f9dSMichael Lange PetscInt numCells = s->last - s->first + 1; 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numCells, (PetscInt**)&s->data)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index==2012 ? PETSC_TRUE : PETSC_FALSE)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(s->data)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 144895fcc3eSMichael Lange } 1452f0bd6dcSMichael Lange } 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1472f0bd6dcSMichael Lange 14819d58f9dSMichael Lange } else if (s->index == 13 || s->index == 2013) { /* Faces */ 1495f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1502f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1512c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1522f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1532f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1542c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 4,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 155f7320561SMichael Lange } else { /* Data section */ 156137d0469SJed Brown PetscInt f, numEntries, numFaces; 1572f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1582c71b3e2SJacob Faibussowitsch PetscCheckFalse(snum != 5,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 160f7320561SMichael Lange switch (s->nd) { 161895fcc3eSMichael Lange case 0: numEntries = PETSC_DETERMINE; break; 162f7320561SMichael Lange case 2: numEntries = 2 + 2; break; /* linear */ 163f7320561SMichael Lange case 3: numEntries = 2 + 3; break; /* triangular */ 164f7320561SMichael Lange case 4: numEntries = 2 + 4; break; /* quadrilateral */ 165f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 166f7320561SMichael Lange } 167f7320561SMichael Lange numFaces = s->last-s->first + 1; 168895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 169895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data)); 171895fcc3eSMichael Lange } 172f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 173895fcc3eSMichael Lange if (s->nd == 0) { 174895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 175137d0469SJed Brown PetscInt numFaceVert = 0; 1765f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE)); 177895fcc3eSMichael Lange if (numEntries == PETSC_DETERMINE) { 178895fcc3eSMichael Lange numEntries = numFaceVert + 2; 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data)); 180895fcc3eSMichael Lange } else { 1812c71b3e2SJacob Faibussowitsch PetscCheckFalse(numEntries != numFaceVert + 2,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 182895fcc3eSMichael Lange } 183895fcc3eSMichael Lange } 1845f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt*)s->data)[f*numEntries]), numEntries, PETSC_INT, s->index==2013 ? PETSC_TRUE : PETSC_FALSE)); 185f7320561SMichael Lange } 186895fcc3eSMichael Lange s->nd = numEntries - 2; 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1882f0bd6dcSMichael Lange } 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1902f0bd6dcSMichael Lange 1912f0bd6dcSMichael Lange } else { /* Unknown section type */ 192060da220SMatthew G. Knepley PetscInt depth = 1; 1932f0bd6dcSMichael Lange do { 1942f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 1955f80ce2aSJacob Faibussowitsch do CHKERRQ(PetscViewerRead(viewer, &(buffer[0]), 1, NULL, PETSC_CHAR)); 1962f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 1972f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 1982f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 1992f0bd6dcSMichael Lange } while (depth > 0); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 2012f0bd6dcSMichael Lange } 2022f0bd6dcSMichael Lange PetscFunctionReturn(0); 2032f0bd6dcSMichael Lange } 2042f0bd6dcSMichael Lange 2052f0bd6dcSMichael Lange /*@C 2062f0bd6dcSMichael Lange DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 2072f0bd6dcSMichael Lange 208d083f849SBarry Smith Collective 2092f0bd6dcSMichael Lange 2102f0bd6dcSMichael Lange Input Parameters: 2112f0bd6dcSMichael Lange + comm - The MPI communicator 2122f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file 2132f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 2142f0bd6dcSMichael Lange 2152f0bd6dcSMichael Lange Output Parameter: 2162f0bd6dcSMichael Lange . dm - The DM object representing the mesh 2172f0bd6dcSMichael Lange 2182f0bd6dcSMichael Lange Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 2192f0bd6dcSMichael Lange 2202f0bd6dcSMichael Lange Level: beginner 2212f0bd6dcSMichael Lange 2222f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate() 2232f0bd6dcSMichael Lange @*/ 2242f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 2252f0bd6dcSMichael Lange { 2262f0bd6dcSMichael Lange PetscMPIInt rank; 227c2fc583bSSatish Balay PetscInt c, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE; 228cf554e2cSMatthew G. Knepley PetscInt numFaces = PETSC_DETERMINE, f, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE; 2293ed8799eSMichael Lange PetscInt *faces = NULL, *cellVertices = NULL, *faceZoneIDs = NULL; 2307368db69SLisandro Dalcin DMLabel faceSets = NULL; 2313f6dc66eSMichael Lange PetscInt d, coordSize; 2323f6dc66eSMichael Lange PetscScalar *coords, *coordsIn = NULL; 2333f6dc66eSMichael Lange PetscSection coordSection; 2343f6dc66eSMichael Lange Vec coordinates; 2352f0bd6dcSMichael Lange 2362f0bd6dcSMichael Lange PetscFunctionBegin; 2375f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 2382f0bd6dcSMichael Lange 239dd400576SPatrick Sanan if (rank == 0) { 2402f0bd6dcSMichael Lange FluentSection s; 2412f0bd6dcSMichael Lange do { 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateFluent_ReadSection(viewer, &s)); 2432f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 244f7320561SMichael Lange dim = s.nd; 2452f0bd6dcSMichael Lange 24619d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 247f7320561SMichael Lange if (s.zoneID == 0) numVertices = s.last; 2483f6dc66eSMichael Lange else { 249*28b400f6SJacob Faibussowitsch PetscCheck(!coordsIn,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 2501a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *) s.data; 2513f6dc66eSMichael Lange } 2522f0bd6dcSMichael Lange 25319d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 254f7320561SMichael Lange if (s.zoneID == 0) numCells = s.last; 255f7320561SMichael Lange else { 256f7320561SMichael Lange switch (s.nd) { 257895fcc3eSMichael Lange case 0: numCellVertices = PETSC_DETERMINE; break; 258f7320561SMichael Lange case 1: numCellVertices = 3; break; /* triangular */ 259f7320561SMichael Lange case 2: numCellVertices = 4; break; /* tetrahedral */ 260f7320561SMichael Lange case 3: numCellVertices = 4; break; /* quadrilateral */ 261f7320561SMichael Lange case 4: numCellVertices = 8; break; /* hexahedral */ 262f7320561SMichael Lange case 5: numCellVertices = 5; break; /* pyramid */ 263f7320561SMichael Lange case 6: numCellVertices = 6; break; /* wedge */ 264895fcc3eSMichael Lange default: numCellVertices = PETSC_DETERMINE; 265f7320561SMichael Lange } 266f7320561SMichael Lange } 2672f0bd6dcSMichael Lange 26819d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 269f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 270930bae4bSMatthew G. Knepley numFaces = (PetscInt) (s.last - s.first + 1); 27135462f7fSMichael Lange if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 272895fcc3eSMichael Lange else numFaceVertices = s.nd; 273f7320561SMichael Lange } else { /* Data section */ 274cf554e2cSMatthew G. Knepley unsigned int z; 275cf554e2cSMatthew G. Knepley 2762c71b3e2SJacob Faibussowitsch PetscCheckFalse(numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 2772c71b3e2SJacob Faibussowitsch PetscCheckFalse(numFaces < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 278f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 279f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 2805f80ce2aSJacob Faibussowitsch if (!faces) CHKERRQ(PetscMalloc1(numFaces*numFaceEntries, &faces)); 2815f80ce2aSJacob Faibussowitsch if (!faceZoneIDs) CHKERRQ(PetscMalloc1(numFaces, &faceZoneIDs)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMemcpy(&faces[(s.first-1)*numFaceEntries], s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt))); 283ec78a56aSMichael Lange /* Record the zoneID for each face set */ 284cf554e2cSMatthew G. Knepley for (z = s.first -1; z < s.last; z++) faceZoneIDs[z] = s.zoneID; 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(s.data)); 286f7320561SMichael Lange } 2872f0bd6dcSMichael Lange } 2882f0bd6dcSMichael Lange } while (s.index >= 0); 2892f0bd6dcSMichael Lange } 2905f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm)); 2912c71b3e2SJacob Faibussowitsch PetscCheckFalse(dim < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 292f7320561SMichael Lange 293f7320561SMichael Lange /* Allocate cell-vertex mesh */ 2945f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetChart(*dm, 0, numCells + numVertices)); 298dd400576SPatrick Sanan if (rank == 0) { 2992c71b3e2SJacob Faibussowitsch PetscCheckFalse(numCells < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 300895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 301895fcc3eSMichael Lange if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1; 3025f80ce2aSJacob Faibussowitsch for (c = 0; c < numCells; ++c) CHKERRQ(DMPlexSetConeSize(*dm, c, numCellVertices)); 303f7320561SMichael Lange } 3045f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(*dm)); 305f7320561SMichael Lange 306dd400576SPatrick Sanan if (rank == 0 && faces) { 307f7320561SMichael Lange /* Derive cell-vertex list from face-vertex and face-cell maps */ 3085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numCells*numCellVertices, &cellVertices)); 309f7320561SMichael Lange for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 310f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 311f7320561SMichael Lange PetscInt *cell; 312f7320561SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 313f7320561SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 314f7320561SMichael Lange const PetscInt *face = &(faces[f*numFaceEntries]); 315f7320561SMichael Lange 316f7320561SMichael Lange if (cl > 0) { 317f7320561SMichael Lange cell = &(cellVertices[(cl-1) * numCellVertices]); 318f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 319f7320561SMichael Lange PetscBool found = PETSC_FALSE; 320f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 321f7320561SMichael Lange if (cell[c] < 0) break; 322f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 323f7320561SMichael Lange } 324f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 325f7320561SMichael Lange } 326f7320561SMichael Lange } 327f7320561SMichael Lange if (cr > 0) { 328f7320561SMichael Lange cell = &(cellVertices[(cr-1) * numCellVertices]); 329f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 330f7320561SMichael Lange PetscBool found = PETSC_FALSE; 331f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 332f7320561SMichael Lange if (cell[c] < 0) break; 333f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 334f7320561SMichael Lange } 335f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 336f7320561SMichael Lange } 337f7320561SMichael Lange } 338f7320561SMichael Lange } 339f7320561SMichael Lange for (c = 0; c < numCells; c++) { 3405f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]))); 341f7320561SMichael Lange } 342f7320561SMichael Lange } 3435f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSymmetrize(*dm)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexStratify(*dm)); 345f7320561SMichael Lange if (interpolate) { 3465fd9971aSMatthew G. Knepley DM idm; 347f7320561SMichael Lange 3485f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 3495f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 350f7320561SMichael Lange *dm = idm; 351f7320561SMichael Lange } 352f7320561SMichael Lange 353dd400576SPatrick Sanan if (rank == 0 && faces) { 354631eb916SMichael Lange PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 355631eb916SMichael Lange const PetscInt *join, *meet; 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(numFaceVertices, &fverts)); 357ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 358ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 359631eb916SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices] - 1; 360631eb916SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1] - 1; 361631eb916SMichael Lange if (cl > 0 && cr > 0) { 362631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 363631eb916SMichael Lange cells[0] = cl; cells[1] = cr; 3645f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet)); 3652c71b3e2SJacob Faibussowitsch PetscCheckFalse(meetSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 3665f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], faceZoneIDs[f])); 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet)); 368631eb916SMichael Lange } else { 369ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1; 3705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 3712c71b3e2SJacob Faibussowitsch PetscCheckFalse(joinSize != 1,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], faceZoneIDs[f])); 3735f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 374ec78a56aSMichael Lange } 375631eb916SMichael Lange } 3765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(fverts)); 377ec78a56aSMichael Lange } 378ec78a56aSMichael Lange 3797368db69SLisandro Dalcin { /* Create Face Sets label at all processes */ 3807368db69SLisandro Dalcin enum {n = 1}; 3817368db69SLisandro Dalcin PetscBool flag[n]; 3827368db69SLisandro Dalcin 3837368db69SLisandro Dalcin flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE; 3845f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 3855f80ce2aSJacob Faibussowitsch if (flag[0]) CHKERRQ(DMCreateLabel(*dm, "Face Sets")); 3867368db69SLisandro Dalcin } 3877368db69SLisandro Dalcin 3883f6dc66eSMichael Lange /* Read coordinates */ 3895f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateSection(*dm, &coordSection)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(coordSection, 1)); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(coordSection, 0, dim)); 3925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 3933f6dc66eSMichael Lange for (v = numCells; v < numCells+numVertices; ++v) { 3945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(coordSection, v, dim)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 3963f6dc66eSMichael Lange } 3975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(coordSection)); 3985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionGetStorageSize(coordSection, &coordSize)); 3995f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_SELF, &coordinates)); 4005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetType(coordinates, VECSTANDARD)); 4035f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 404dd400576SPatrick Sanan if (rank == 0 && coordsIn) { 4053f6dc66eSMichael Lange for (v = 0; v < numVertices; ++v) { 4063f6dc66eSMichael Lange for (d = 0; d < dim; ++d) { 4073f6dc66eSMichael Lange coords[v*dim+d] = coordsIn[v*dim+d]; 4083f6dc66eSMichael Lange } 4093f6dc66eSMichael Lange } 4103f6dc66eSMichael Lange } 4115f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 4125f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinatesLocal(*dm, coordinates)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&coordinates)); 4147368db69SLisandro Dalcin 415dd400576SPatrick Sanan if (rank == 0) { 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cellVertices)); 4175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(faces)); 4185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(faceZoneIDs)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(coordsIn)); 420f7320561SMichael Lange } 4212f0bd6dcSMichael Lange PetscFunctionReturn(0); 4222f0bd6dcSMichael Lange } 423