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 55733454dSMatthew G. Knepley /* Utility struct to store the contents of a Fluent file in memory */ 65733454dSMatthew G. Knepley typedef struct { 75733454dSMatthew G. Knepley int index; /* Type of section */ 85733454dSMatthew G. Knepley unsigned int zoneID; 95733454dSMatthew G. Knepley unsigned int first; 105733454dSMatthew G. Knepley unsigned int last; 115733454dSMatthew G. Knepley int type; 125733454dSMatthew G. Knepley int nd; /* Either ND or element-type */ 135733454dSMatthew G. Knepley void *data; 145733454dSMatthew G. Knepley } FluentSection; 155733454dSMatthew G. Knepley 16cc4c1da9SBarry Smith /*@ 17a1cb98faSBarry Smith DMPlexCreateFluentFromFile - Create a `DMPLEX` mesh from a Fluent mesh file 182f0bd6dcSMichael Lange 19a1cb98faSBarry Smith Collective 20a1cb98faSBarry Smith 21a1cb98faSBarry Smith Input Parameters: 222f0bd6dcSMichael Lange + comm - The MPI communicator 232f0bd6dcSMichael Lange . filename - Name of the Fluent mesh file 242f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 252f0bd6dcSMichael Lange 262f0bd6dcSMichael Lange Output Parameter: 27a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 282f0bd6dcSMichael Lange 292f0bd6dcSMichael Lange Level: beginner 302f0bd6dcSMichael Lange 311cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateFluent()`, `DMPlexCreate()` 322f0bd6dcSMichael Lange @*/ 33d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 34d71ae5a4SJacob Faibussowitsch { 352f0bd6dcSMichael Lange PetscViewer viewer; 362f0bd6dcSMichael Lange 372f0bd6dcSMichael Lange PetscFunctionBegin; 382f0bd6dcSMichael Lange /* Create file viewer and build plex */ 399566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 409566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 419566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 439566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent(comm, viewer, interpolate, dm)); 449566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 462f0bd6dcSMichael Lange } 472f0bd6dcSMichael Lange 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 49d71ae5a4SJacob Faibussowitsch { 504b375a39SMichael Lange PetscInt ret, i = 0; 512f0bd6dcSMichael Lange 522f0bd6dcSMichael Lange PetscFunctionBegin; 53f4f49eeaSPierre Jolivet do PetscCall(PetscViewerRead(viewer, &buffer[i++], 1, &ret, PETSC_CHAR)); 54c7cc6ee6SMatthew G. Knepley while (ret > 0 && buffer[i - 1] != '\0' && buffer[i - 1] != delim && i < PETSC_MAX_PATH_LEN - 1); 559371c9d4SSatish Balay if (!ret) buffer[i - 1] = '\0'; 569371c9d4SSatish Balay else buffer[i] = '\0'; 57c7cc6ee6SMatthew G. Knepley PetscCheck(i < PETSC_MAX_PATH_LEN - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Buffer overflow! This is not a valid Fluent file."); 583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 592f0bd6dcSMichael Lange } 602f0bd6dcSMichael Lange 61c7cc6ee6SMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary, PetscInt *numClosingParens) 62d71ae5a4SJacob Faibussowitsch { 6370c9a859SSatish Balay int fdes = 0; 6419d58f9dSMichael Lange FILE *file; 6519d58f9dSMichael Lange PetscInt i; 6619d58f9dSMichael Lange 6719d58f9dSMichael Lange PetscFunctionBegin; 68c7cc6ee6SMatthew G. Knepley *numClosingParens = 0; 6919d58f9dSMichael Lange if (binary) { 7019d58f9dSMichael Lange /* Extract raw file descriptor to read binary block */ 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetPointer(viewer, &file)); 72c69effb2SJacob Faibussowitsch PetscCall(PetscFFlush(file)); 739371c9d4SSatish Balay fdes = fileno(file); 7419d58f9dSMichael Lange } 7519d58f9dSMichael Lange 7619d58f9dSMichael Lange if (!binary && dtype == PETSC_INT) { 7719d58f9dSMichael Lange char cbuf[256]; 78cfb60857SMatthew G. Knepley unsigned int ibuf; 79cfb60857SMatthew G. Knepley int snum; 8019d58f9dSMichael Lange /* Parse hexadecimal ascii integers */ 8119d58f9dSMichael Lange for (i = 0; i < count; i++) { 82c7cc6ee6SMatthew G. Knepley size_t len; 83c7cc6ee6SMatthew G. Knepley 849566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING)); 8519d58f9dSMichael Lange snum = sscanf(cbuf, "%x", &ibuf); 8608401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 8719d58f9dSMichael Lange ((PetscInt *)data)[i] = (PetscInt)ibuf; 88c7cc6ee6SMatthew G. Knepley // Check for trailing parentheses 89c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrlen(cbuf, &len)); 90c7cc6ee6SMatthew G. Knepley while (cbuf[len - 1] == ')' && len > 0) { 91c7cc6ee6SMatthew G. Knepley ++(*numClosingParens); 92c7cc6ee6SMatthew G. Knepley --len; 93c7cc6ee6SMatthew G. Knepley } 9419d58f9dSMichael Lange } 9519d58f9dSMichael Lange } else if (binary && dtype == PETSC_INT) { 9619d58f9dSMichael Lange /* Always read 32-bit ints and cast to PetscInt */ 9719d58f9dSMichael Lange int *ibuf; 989566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, &ibuf)); 999566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM)); 1009566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, count)); 101835f2295SStefano Zampini for (i = 0; i < count; i++) ((PetscInt *)data)[i] = ibuf[i]; 1029566063dSJacob Faibussowitsch PetscCall(PetscFree(ibuf)); 10319d58f9dSMichael Lange 10419d58f9dSMichael Lange } else if (binary && dtype == PETSC_SCALAR) { 10519d58f9dSMichael Lange float *fbuf; 10619d58f9dSMichael Lange /* Always read 32-bit floats and cast to PetscScalar */ 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, &fbuf)); 1089566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT)); 1099566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(fbuf, PETSC_FLOAT, count)); 110835f2295SStefano Zampini for (i = 0; i < count; i++) ((PetscScalar *)data)[i] = fbuf[i]; 1119566063dSJacob Faibussowitsch PetscCall(PetscFree(fbuf)); 11219d58f9dSMichael Lange } else { 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIRead(viewer, data, count, NULL, dtype)); 11419d58f9dSMichael Lange } 1153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11619d58f9dSMichael Lange } 11719d58f9dSMichael Lange 118d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 119d71ae5a4SJacob Faibussowitsch { 1202f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 1212f0bd6dcSMichael Lange int snum; 1222f0bd6dcSMichael Lange 1232f0bd6dcSMichael Lange PetscFunctionBegin; 1242f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 1259566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1269566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ' ')); 127f4f49eeaSPierre Jolivet snum = sscanf(buffer, "%d", &s->index); 1282f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 1299371c9d4SSatish Balay if (snum < 1) { 1309371c9d4SSatish Balay s->index = -1; 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1329371c9d4SSatish Balay } 1332f0bd6dcSMichael Lange 1342f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 1359566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1362f0bd6dcSMichael Lange 1372f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 1389566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 139f4f49eeaSPierre Jolivet snum = sscanf(buffer, "%d", &s->nd); 14008401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1412f0bd6dcSMichael Lange 14219d58f9dSMichael Lange } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 143c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 144c7cc6ee6SMatthew G. Knepley 1459566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 146f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 14708401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1483f6dc66eSMichael Lange if (s->zoneID > 0) { 1493f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 1509566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s->nd * numCoords, (PetscScalar **)&s->data)); 152c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd * numCoords, PETSC_SCALAR, s->index == 2010 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 153c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 154c7cc6ee6SMatthew G. Knepley else --numClosingParens; 1553f6dc66eSMichael Lange } 156c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 157c7cc6ee6SMatthew G. Knepley else --numClosingParens; 158c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 15919d58f9dSMichael Lange } else if (s->index == 12 || s->index == 2012) { /* Cells */ 160c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 161c7cc6ee6SMatthew G. Knepley 1629566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 163f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x", &s->zoneID); 16408401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1652f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 166f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd); 16708401ef6SPierre Jolivet PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 168f7320561SMichael Lange } else { /* Data section */ 169f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 17008401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 171895fcc3eSMichael Lange if (s->nd == 0) { 172895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 17319d58f9dSMichael Lange PetscInt numCells = s->last - s->first + 1; 1749566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, (PetscInt **)&s->data)); 176c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index == 2012 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 177c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 178c7cc6ee6SMatthew G. Knepley else --numClosingParens; 179895fcc3eSMichael Lange } 1802f0bd6dcSMichael Lange } 181c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 182c7cc6ee6SMatthew G. Knepley else --numClosingParens; 183c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 18419d58f9dSMichael Lange } else if (s->index == 13 || s->index == 2013) { /* Faces */ 185c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 186c7cc6ee6SMatthew G. Knepley 1879566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 188f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x", &s->zoneID); 18908401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1902f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 191f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd); 19208401ef6SPierre Jolivet PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 193f7320561SMichael Lange } else { /* Data section */ 1945733454dSMatthew G. Knepley PetscInt numEntries, numFaces, maxsize = 0, offset = 0; 1955733454dSMatthew G. Knepley 196f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 19708401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1989566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 199f7320561SMichael Lange switch (s->nd) { 200d71ae5a4SJacob Faibussowitsch case 0: 201d71ae5a4SJacob Faibussowitsch numEntries = PETSC_DETERMINE; 202d71ae5a4SJacob Faibussowitsch break; 203d71ae5a4SJacob Faibussowitsch case 2: 204d71ae5a4SJacob Faibussowitsch numEntries = 2 + 2; 205d71ae5a4SJacob Faibussowitsch break; /* linear */ 206d71ae5a4SJacob Faibussowitsch case 3: 207d71ae5a4SJacob Faibussowitsch numEntries = 2 + 3; 208d71ae5a4SJacob Faibussowitsch break; /* triangular */ 209d71ae5a4SJacob Faibussowitsch case 4: 210d71ae5a4SJacob Faibussowitsch numEntries = 2 + 4; 211d71ae5a4SJacob Faibussowitsch break; /* quadrilateral */ 212d71ae5a4SJacob Faibussowitsch default: 213d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 214f7320561SMichael Lange } 215f7320561SMichael Lange numFaces = s->last - s->first + 1; 216895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 217895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 2189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data)); 219895fcc3eSMichael Lange } 2205733454dSMatthew G. Knepley for (PetscInt f = 0; f < numFaces; f++) { 221895fcc3eSMichael Lange if (s->nd == 0) { 222895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 223137d0469SJed Brown PetscInt numFaceVert = 0; 224c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 2255733454dSMatthew G. Knepley if (!f) { 2265733454dSMatthew G. Knepley maxsize = (numFaceVert + 3) * numFaces; 2275733454dSMatthew G. Knepley PetscCall(PetscMalloc1(maxsize, (PetscInt **)&s->data)); 228895fcc3eSMichael Lange } else { 2295733454dSMatthew G. Knepley if (offset + numFaceVert + 3 >= maxsize) { 2305733454dSMatthew G. Knepley PetscInt *tmp; 2315733454dSMatthew G. Knepley 2325733454dSMatthew G. Knepley PetscCall(PetscMalloc1(maxsize * 2, &tmp)); 2335733454dSMatthew G. Knepley PetscCall(PetscArraycpy(tmp, (PetscInt *)s->data, maxsize)); 2345733454dSMatthew G. Knepley PetscCall(PetscFree(s->data)); 2355733454dSMatthew G. Knepley maxsize *= 2; 2365733454dSMatthew G. Knepley s->data = tmp; 237895fcc3eSMichael Lange } 238895fcc3eSMichael Lange } 2395733454dSMatthew G. Knepley ((PetscInt *)s->data)[offset] = numFaceVert; 2405733454dSMatthew G. Knepley ++offset; 2415733454dSMatthew G. Knepley numEntries = numFaceVert + 2; 242f7320561SMichael Lange } 2435733454dSMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt *)s->data)[offset]), numEntries, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 2445733454dSMatthew G. Knepley offset += numEntries; 2455733454dSMatthew G. Knepley } 2465733454dSMatthew G. Knepley if (s->nd != 0) PetscCall(PetscMPIIntCast(numEntries - 2, &s->nd)); 247c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 248c7cc6ee6SMatthew G. Knepley else --numClosingParens; 2492f0bd6dcSMichael Lange } 250c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 251c7cc6ee6SMatthew G. Knepley else --numClosingParens; 252c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 253c7cc6ee6SMatthew G. Knepley } else if (s->index == 39) { /* Label information */ 254c7cc6ee6SMatthew G. Knepley char labelName[PETSC_MAX_PATH_LEN]; 255c7cc6ee6SMatthew G. Knepley char caseName[PETSC_MAX_PATH_LEN]; 2562f0bd6dcSMichael Lange 257c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 258c3279546SPierre Jolivet snum = sscanf(buffer, "(%u %s %s %d)", &s->zoneID, caseName, labelName, &s->nd); 259c7cc6ee6SMatthew G. Knepley PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file: %d", snum); 260c7cc6ee6SMatthew G. Knepley PetscInt depth = 1; 261c7cc6ee6SMatthew G. Knepley do { 262c7cc6ee6SMatthew G. Knepley /* Match parentheses when parsing unknown sections */ 263c7cc6ee6SMatthew G. Knepley do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR)); 264c7cc6ee6SMatthew G. Knepley while (buffer[0] != '(' && buffer[0] != ')'); 265c7cc6ee6SMatthew G. Knepley if (buffer[0] == '(') depth++; 266c7cc6ee6SMatthew G. Knepley if (buffer[0] == ')') depth--; 267c7cc6ee6SMatthew G. Knepley } while (depth > 0); 268c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 269c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrallocpy(labelName, (char **)&s->data)); 270c3279546SPierre Jolivet PetscCall(PetscInfo((PetscObject)viewer, "CASE: Zone ID %u is label %s\n", s->zoneID, labelName)); 2712f0bd6dcSMichael Lange } else { /* Unknown section type */ 272060da220SMatthew G. Knepley PetscInt depth = 1; 2732f0bd6dcSMichael Lange do { 2742f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 275f4f49eeaSPierre Jolivet do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR)); 2762f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 2772f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 2782f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 2792f0bd6dcSMichael Lange } while (depth > 0); 2809566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 2812f0bd6dcSMichael Lange } 2823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2832f0bd6dcSMichael Lange } 2842f0bd6dcSMichael Lange 285c7cc6ee6SMatthew G. Knepley // Inserts point `face` with orientation `ornt` into the cone of point `cell` at position `c`, which is the first empty slot 286c7cc6ee6SMatthew G. Knepley static PetscErrorCode InsertFace(DM dm, PetscInt cell, PetscInt face, PetscInt ornt) 287c7cc6ee6SMatthew G. Knepley { 288c7cc6ee6SMatthew G. Knepley const PetscInt *cone; 289c7cc6ee6SMatthew G. Knepley PetscInt coneSize, c; 290c7cc6ee6SMatthew G. Knepley 291c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 292c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cell, &cone)); 293c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 294c7cc6ee6SMatthew G. Knepley for (c = 0; c < coneSize; ++c) 295c7cc6ee6SMatthew G. Knepley if (cone[c] < 0) break; 2965733454dSMatthew G. Knepley PetscCheck(c < coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Face %" PetscInt_FMT " could not be inserted in cone of cell %" PetscInt_FMT " with size %" PetscInt_FMT, face, cell, coneSize); 297c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertCone(dm, cell, c, face)); 298c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertConeOrientation(dm, cell, c, ornt)); 299c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 300c7cc6ee6SMatthew G. Knepley } 301c7cc6ee6SMatthew G. Knepley 302c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderPolygon(DM dm, PetscInt cell) 303c7cc6ee6SMatthew G. Knepley { 304c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt; 305c7cc6ee6SMatthew G. Knepley PetscInt coneSize, newCone[16], newOrnt[16]; 306c7cc6ee6SMatthew G. Knepley 307c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 308c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 309c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 310c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 311c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 312c7cc6ee6SMatthew G. Knepley for (PetscInt c = 1; c < coneSize; ++c) { 313c7cc6ee6SMatthew G. Knepley const PetscInt *fcone; 314c7cc6ee6SMatthew G. Knepley PetscInt firstVertex, lastVertex, c2; 315c7cc6ee6SMatthew G. Knepley 316c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[c - 1], &fcone)); 317c7cc6ee6SMatthew G. Knepley lastVertex = newOrnt[c - 1] ? fcone[0] : fcone[1]; 318c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < coneSize; ++c2) { 319c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2; 320c7cc6ee6SMatthew G. Knepley 321c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cone[c2], &fcone2)); 322c7cc6ee6SMatthew G. Knepley firstVertex = ornt[c2] ? fcone2[1] : fcone2[0]; 323c7cc6ee6SMatthew G. Knepley if (lastVertex == firstVertex) { 324c7cc6ee6SMatthew G. Knepley // Point `cell` matched point `lastVertex` on face `cone[c2]` with orientation `ornt[c2]` 325c7cc6ee6SMatthew G. Knepley break; 326c7cc6ee6SMatthew G. Knepley } 327c7cc6ee6SMatthew G. Knepley } 328c7cc6ee6SMatthew G. Knepley PetscCheck(c2 < coneSize, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match as position %" PetscInt_FMT, cell, c); 329c7cc6ee6SMatthew G. Knepley newCone[c] = cone[c2]; 330c7cc6ee6SMatthew G. Knepley newOrnt[c] = ornt[c2]; 331c7cc6ee6SMatthew G. Knepley } 332c7cc6ee6SMatthew G. Knepley { 333c7cc6ee6SMatthew G. Knepley const PetscInt *fcone, *fcone2; 334c7cc6ee6SMatthew G. Knepley PetscInt vertex, vertex2; 335c7cc6ee6SMatthew G. Knepley 336c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[coneSize - 1], &fcone)); 337c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[0], &fcone2)); 338c7cc6ee6SMatthew G. Knepley vertex = newOrnt[coneSize - 1] ? fcone[0] : fcone[1]; 339c7cc6ee6SMatthew G. Knepley vertex2 = newOrnt[0] ? fcone2[1] : fcone2[0]; 340c7cc6ee6SMatthew G. Knepley PetscCheck(vertex == vertex2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " did not match at the endpoint", cell); 341c7cc6ee6SMatthew G. Knepley } 342c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 343c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 344c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 345c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 346c7cc6ee6SMatthew G. Knepley } 347c7cc6ee6SMatthew G. Knepley 348fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderTetrahedron(PetscViewer viewer, DM dm, PetscInt cell) 349c7cc6ee6SMatthew G. Knepley { 350c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr, faces[4] = {0, 1, 3, 2}; 351c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 352c7cc6ee6SMatthew G. Knepley 353c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 354c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 355c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 356c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 357c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 358c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]); 359c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial triangle 360c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 3; ++e) { 361c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 362c7cc6ee6SMatthew G. Knepley PetscInt c; 363c7cc6ee6SMatthew G. Knepley 364c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the tetrahedron 365c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 366c7cc6ee6SMatthew G. Knepley for (c = 1; c < 4; ++c) { 367c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 368c7cc6ee6SMatthew G. Knepley PetscInt c2; 369fc8f8d65SMatthew G. Knepley PetscBool flip = PETSC_FALSE; 370c7cc6ee6SMatthew G. Knepley 371c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 372c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 373c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]); 374c7cc6ee6SMatthew G. Knepley // Check for edge 375c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 3; ++c2) { 376fc8f8d65SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 377c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 378c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 379fc8f8d65SMatthew G. Knepley //PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell % " PetscInt_FMT " edge %" PetscInt_FMT " (%" PetscInt_FMT ") found twice with the same orientation in face %" PetscInt_FMT " edge %" PetscInt_FMT, cell, edge, e, c, c2); 380c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 381fc8f8d65SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Matched cell %" PetscInt_FMT " edge %" PetscInt_FMT "/%" PetscInt_FMT " (%" PetscInt_FMT ") to face %" PetscInt_FMT "/%" PetscInt_FMT " edge %" PetscInt_FMT " (%" PetscInt_FMT ")\n", cell, edge, e, eornt, cone[c], c, c2, eornt2)); 382fc8f8d65SMatthew G. Knepley flip = eornt != -(eornt2 + 1) ? PETSC_TRUE : PETSC_FALSE; 383c7cc6ee6SMatthew G. Knepley break; 384c7cc6ee6SMatthew G. Knepley } 385c7cc6ee6SMatthew G. Knepley } 386c7cc6ee6SMatthew G. Knepley if (c2 < 3) { 387c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 388c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 389fc8f8d65SMatthew G. Knepley // Face 1 should match its edge 2 390fc8f8d65SMatthew G. Knepley // Face 2 should match its edge 0 391fc8f8d65SMatthew G. Knepley // Face 3 should match its edge 0 392fc8f8d65SMatthew G. Knepley if (flip) { 393fc8f8d65SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, -((c2 + (!e ? 1 : 2)) % 3 + 1), ornt[c]); 394fc8f8d65SMatthew G. Knepley } else { 395c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, !e ? (c2 + 1) % 3 : c2, ornt[c]); 396fc8f8d65SMatthew G. Knepley } 397c7cc6ee6SMatthew G. Knepley break; 398c7cc6ee6SMatthew G. Knepley } 399c7cc6ee6SMatthew G. Knepley } 400c7cc6ee6SMatthew G. Knepley PetscCheck(c < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e); 401c7cc6ee6SMatthew G. Knepley } 402c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 403c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 404c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 405c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 406c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 407c7cc6ee6SMatthew G. Knepley } 408c7cc6ee6SMatthew G. Knepley 409c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderHexahedron(DM dm, PetscInt cell) 410c7cc6ee6SMatthew G. Knepley { 411c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr; 412c7cc6ee6SMatthew G. Knepley const PetscInt faces[6] = {0, 5, 3, 4, 2, 1}; 413c7cc6ee6SMatthew G. Knepley PetscInt used[6] = {1, 0, 0, 0, 0, 0}; 414c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 415c7cc6ee6SMatthew G. Knepley 416c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 417c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 418c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 419c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 420c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 421c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[0]); 422c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial quadrilateral 423c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 4; ++e) { 424c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 425c7cc6ee6SMatthew G. Knepley PetscInt c; 426c7cc6ee6SMatthew G. Knepley 427c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the hexahedron 428c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 429c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) { 430c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 431c7cc6ee6SMatthew G. Knepley PetscInt c2; 432c7cc6ee6SMatthew G. Knepley 433c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 434c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 435c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 436c7cc6ee6SMatthew G. Knepley // Check for edge 437c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 438c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 439c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 440c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 441c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 442c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 443c7cc6ee6SMatthew G. Knepley break; 444c7cc6ee6SMatthew G. Knepley } 445c7cc6ee6SMatthew G. Knepley } 446c7cc6ee6SMatthew G. Knepley if (c2 < 4) { 447c7cc6ee6SMatthew G. Knepley used[c] = 1; 448c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 449c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 450c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, !e ? (c2 + 1) % 4 : c2, ornt[c]); 451c7cc6ee6SMatthew G. Knepley break; 452c7cc6ee6SMatthew G. Knepley } 453c7cc6ee6SMatthew G. Knepley } 454c7cc6ee6SMatthew G. Knepley PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e); 455c7cc6ee6SMatthew G. Knepley } 456c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 457c7cc6ee6SMatthew G. Knepley // Add last face 458c7cc6ee6SMatthew G. Knepley { 459c7cc6ee6SMatthew G. Knepley PetscInt c, c2; 460c7cc6ee6SMatthew G. Knepley 461c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) 462c7cc6ee6SMatthew G. Knepley if (!used[c]) break; 463c7cc6ee6SMatthew G. Knepley PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell); 464c7cc6ee6SMatthew G. Knepley // Match first edge to 3rd edge in newCone[2] 465c7cc6ee6SMatthew G. Knepley { 466c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 467c7cc6ee6SMatthew G. Knepley 468c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt)); 469c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]); 470c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 471c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 472c7cc6ee6SMatthew G. Knepley 473c7cc6ee6SMatthew G. Knepley const PetscInt e = 2; 474c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 475c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge` with final orientation `eornt` of face `newCone[2]` to some edge of face `cone[c]` with orientation `ornt[c]` 476c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 477c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 478c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 479c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 480c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 481c7cc6ee6SMatthew G. Knepley // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 482c7cc6ee6SMatthew G. Knepley break; 483c7cc6ee6SMatthew G. Knepley } 484c7cc6ee6SMatthew G. Knepley } 485c7cc6ee6SMatthew G. Knepley PetscCheck(c2 < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in"); 486c7cc6ee6SMatthew G. Knepley } 487c7cc6ee6SMatthew G. Knepley newCone[faces[5]] = cone[c]; 488c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched 489c7cc6ee6SMatthew G. Knepley newOrnt[faces[5]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]); 490c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 491c7cc6ee6SMatthew G. Knepley } 492c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 493c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 494c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 495c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 496c7cc6ee6SMatthew G. Knepley } 497c7cc6ee6SMatthew G. Knepley 4985733454dSMatthew G. Knepley // {0, 1, 2}, {3, 4, 5}, {0, 2, 4, 3}, {2, 1, 5, 4}, {1, 0, 3, 5} 4995733454dSMatthew G. Knepley static PetscErrorCode ReorderWedge(DM dm, PetscInt cell) 5005733454dSMatthew G. Knepley { 5015733454dSMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr; 5025733454dSMatthew G. Knepley const PetscInt faces[5] = {0, 4, 3, 2, 1}; 5035733454dSMatthew G. Knepley PetscInt used[5] = {0, 0, 0, 0, 0}; 5045733454dSMatthew G. Knepley PetscInt newCone[16], newOrnt[16], cS, bottom = 0; 5055733454dSMatthew G. Knepley 5065733454dSMatthew G. Knepley PetscFunctionBegin; 5075733454dSMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &cS)); 5085733454dSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 5095733454dSMatthew G. Knepley for (PetscInt c = 0; c < cS; ++c) { 5105733454dSMatthew G. Knepley DMPolytopeType ct; 5115733454dSMatthew G. Knepley 5125733454dSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cone[c], &ct)); 5135733454dSMatthew G. Knepley if (ct == DM_POLYTOPE_TRIANGLE) { 5145733454dSMatthew G. Knepley bottom = c; 5155733454dSMatthew G. Knepley break; 5165733454dSMatthew G. Knepley } 5175733454dSMatthew G. Knepley } 5185733454dSMatthew G. Knepley used[bottom] = 1; 5195733454dSMatthew G. Knepley newCone[0] = cone[bottom]; 5205733454dSMatthew G. Knepley newOrnt[0] = ornt[bottom]; 5215733454dSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 5225733454dSMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]); 5235733454dSMatthew G. Knepley // Loop over each edge in the initial triangle 5245733454dSMatthew G. Knepley for (PetscInt e = 0; e < 3; ++e) { 5255733454dSMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 5265733454dSMatthew G. Knepley PetscInt c; 5275733454dSMatthew G. Knepley 5285733454dSMatthew G. Knepley // Loop over each remaining face in the prism 5295733454dSMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 5305733454dSMatthew G. Knepley for (c = 0; c < 5; ++c) { 5315733454dSMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 5325733454dSMatthew G. Knepley DMPolytopeType ct; 5335733454dSMatthew G. Knepley PetscInt c2; 5345733454dSMatthew G. Knepley 5355733454dSMatthew G. Knepley if (c == bottom) continue; 5365733454dSMatthew G. Knepley PetscCall(DMPlexGetCellType(dm, cone[c], &ct)); 5375733454dSMatthew G. Knepley if (ct != DM_POLYTOPE_QUADRILATERAL) continue; 5385733454dSMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 5395733454dSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 5405733454dSMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 5415733454dSMatthew G. Knepley // Check for edge 5425733454dSMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 5435733454dSMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 5445733454dSMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 5455733454dSMatthew G. Knepley if (edge == edge2) { 5465733454dSMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 5475733454dSMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 5485733454dSMatthew G. Knepley break; 5495733454dSMatthew G. Knepley } 5505733454dSMatthew G. Knepley } 5515733454dSMatthew G. Knepley if (c2 < 4) { 5525733454dSMatthew G. Knepley used[c] = 1; 5535733454dSMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 5545733454dSMatthew G. Knepley // Compute new orientation of face based on which edge was matched, edge 0 should always match the bottom 5555733454dSMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]); 5565733454dSMatthew G. Knepley break; 5575733454dSMatthew G. Knepley } 5585733454dSMatthew G. Knepley } 5595733454dSMatthew G. Knepley PetscCheck(c < 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find a face match for edge %" PetscInt_FMT, cell, e); 5605733454dSMatthew G. Knepley } 5615733454dSMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 5625733454dSMatthew G. Knepley // Add last face 5635733454dSMatthew G. Knepley { 5645733454dSMatthew G. Knepley PetscInt c, c2; 5655733454dSMatthew G. Knepley 5665733454dSMatthew G. Knepley for (c = 0; c < 5; ++c) 5675733454dSMatthew G. Knepley if (!used[c]) break; 5685733454dSMatthew G. Knepley PetscCheck(c < 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell); 5695733454dSMatthew G. Knepley // Match first edge to 3rd edge in newCone[2] 5705733454dSMatthew G. Knepley { 5715733454dSMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 5725733454dSMatthew G. Knepley 5735733454dSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt)); 5745733454dSMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]); 5755733454dSMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 5765733454dSMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]); 5775733454dSMatthew G. Knepley 5785733454dSMatthew G. Knepley const PetscInt e = 2; 5795733454dSMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 5805733454dSMatthew G. Knepley // Trying to match edge `edge` with final orientation `eornt` of face `newCone[2]` to some edge of face `cone[c]` with orientation `ornt[c]` 5815733454dSMatthew G. Knepley for (c2 = 0; c2 < 3; ++c2) { 5825733454dSMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 5835733454dSMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 5845733454dSMatthew G. Knepley if (edge == edge2) { 5855733454dSMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 5865733454dSMatthew G. Knepley // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 5875733454dSMatthew G. Knepley break; 5885733454dSMatthew G. Knepley } 5895733454dSMatthew G. Knepley } 5905733454dSMatthew G. Knepley PetscCheck(c2 < 3, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in"); 5915733454dSMatthew G. Knepley } 5925733454dSMatthew G. Knepley newCone[faces[4]] = cone[c]; 5935733454dSMatthew G. Knepley // Compute new orientation of face based on which edge was matched 5945733454dSMatthew G. Knepley newOrnt[faces[4]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, c2, ornt[c]); 5955733454dSMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 5965733454dSMatthew G. Knepley } 5975733454dSMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 5985733454dSMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 5995733454dSMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 6005733454dSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 6015733454dSMatthew G. Knepley } 6025733454dSMatthew G. Knepley 603fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderCell(PetscViewer viewer, DM dm, PetscInt cell, DMPolytopeType ct) 604c7cc6ee6SMatthew G. Knepley { 605c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 606c7cc6ee6SMatthew G. Knepley switch (ct) { 607c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 608c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 609c7cc6ee6SMatthew G. Knepley PetscCall(ReorderPolygon(dm, cell)); 610c7cc6ee6SMatthew G. Knepley break; 611c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TETRAHEDRON: 612fc8f8d65SMatthew G. Knepley PetscCall(ReorderTetrahedron(viewer, dm, cell)); 613c7cc6ee6SMatthew G. Knepley break; 614c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_HEXAHEDRON: 615c7cc6ee6SMatthew G. Knepley PetscCall(ReorderHexahedron(dm, cell)); 616c7cc6ee6SMatthew G. Knepley break; 6175733454dSMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 6185733454dSMatthew G. Knepley PetscCall(ReorderWedge(dm, cell)); 6195733454dSMatthew G. Knepley break; 620c7cc6ee6SMatthew G. Knepley default: 621c7cc6ee6SMatthew G. Knepley PetscCheck(0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Celltype %s is unsupported", DMPolytopeTypes[ct]); 622c7cc6ee6SMatthew G. Knepley break; 623c7cc6ee6SMatthew G. Knepley } 624c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 625c7cc6ee6SMatthew G. Knepley } 626c7cc6ee6SMatthew G. Knepley 6275733454dSMatthew G. Knepley static PetscErrorCode GetNumCellFaces(int nd, PetscInt *numCellFaces, DMPolytopeType *ct) 6285733454dSMatthew G. Knepley { 6295733454dSMatthew G. Knepley PetscFunctionBegin; 6305733454dSMatthew G. Knepley *ct = DM_POLYTOPE_POINT; 6315733454dSMatthew G. Knepley switch (nd) { 6325733454dSMatthew G. Knepley case 0: 6335733454dSMatthew G. Knepley *numCellFaces = PETSC_DETERMINE; 6345733454dSMatthew G. Knepley break; 6355733454dSMatthew G. Knepley case 1: 6365733454dSMatthew G. Knepley *numCellFaces = 3; 6375733454dSMatthew G. Knepley *ct = DM_POLYTOPE_TRIANGLE; 6385733454dSMatthew G. Knepley break; 6395733454dSMatthew G. Knepley case 2: 6405733454dSMatthew G. Knepley *numCellFaces = 4; 6415733454dSMatthew G. Knepley *ct = DM_POLYTOPE_TETRAHEDRON; 6425733454dSMatthew G. Knepley break; 6435733454dSMatthew G. Knepley case 3: 6445733454dSMatthew G. Knepley *numCellFaces = 4; 6455733454dSMatthew G. Knepley *ct = DM_POLYTOPE_QUADRILATERAL; 6465733454dSMatthew G. Knepley break; 6475733454dSMatthew G. Knepley case 4: 6485733454dSMatthew G. Knepley *numCellFaces = 6; 6495733454dSMatthew G. Knepley *ct = DM_POLYTOPE_HEXAHEDRON; 6505733454dSMatthew G. Knepley break; 6515733454dSMatthew G. Knepley case 5: 6525733454dSMatthew G. Knepley *numCellFaces = 5; 6535733454dSMatthew G. Knepley *ct = DM_POLYTOPE_PYRAMID; 6545733454dSMatthew G. Knepley break; 6555733454dSMatthew G. Knepley case 6: 6565733454dSMatthew G. Knepley *numCellFaces = 5; 6575733454dSMatthew G. Knepley *ct = DM_POLYTOPE_TRI_PRISM; 6585733454dSMatthew G. Knepley break; 6595733454dSMatthew G. Knepley default: 6605733454dSMatthew G. Knepley *numCellFaces = PETSC_DETERMINE; 6615733454dSMatthew G. Knepley } 6625733454dSMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 6635733454dSMatthew G. Knepley } 6645733454dSMatthew G. Knepley 6652f0bd6dcSMichael Lange /*@C 6661d27aa22SBarry Smith DMPlexCreateFluent - Create a `DMPLEX` mesh from a Fluent mesh file <http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm>. 6672f0bd6dcSMichael Lange 668d083f849SBarry Smith Collective 6692f0bd6dcSMichael Lange 6702f0bd6dcSMichael Lange Input Parameters: 6712f0bd6dcSMichael Lange + comm - The MPI communicator 672a1cb98faSBarry Smith . viewer - The `PetscViewer` associated with a Fluent mesh file 6732f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 6742f0bd6dcSMichael Lange 6752f0bd6dcSMichael Lange Output Parameter: 676a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 6772f0bd6dcSMichael Lange 6782f0bd6dcSMichael Lange Level: beginner 6792f0bd6dcSMichael Lange 6801cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()` 6812f0bd6dcSMichael Lange @*/ 682d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 683d71ae5a4SJacob Faibussowitsch { 684c7cc6ee6SMatthew G. Knepley PetscInt dim = PETSC_DETERMINE; 685c7cc6ee6SMatthew G. Knepley PetscInt numCells = 0; 686c7cc6ee6SMatthew G. Knepley PetscInt numVertices = 0; 6875733454dSMatthew G. Knepley PetscInt *cellSizes = NULL; 6885733454dSMatthew G. Knepley DMPolytopeType *cellTypes = NULL; 689c7cc6ee6SMatthew G. Knepley PetscInt numFaces = 0; 690c7cc6ee6SMatthew G. Knepley PetscInt *faces = NULL; 6915733454dSMatthew G. Knepley PetscInt *faceSizes = NULL; 6925733454dSMatthew G. Knepley PetscInt *faceAdjCell = NULL; 693c7cc6ee6SMatthew G. Knepley PetscInt *cellVertices = NULL; 694c7cc6ee6SMatthew G. Knepley unsigned int *faceZoneIDs = NULL; 6957368db69SLisandro Dalcin DMLabel faceSets = NULL; 696c7cc6ee6SMatthew G. Knepley DMLabel *zoneLabels = NULL; 697c7cc6ee6SMatthew G. Knepley const char **zoneNames = NULL; 698c7cc6ee6SMatthew G. Knepley unsigned int maxZoneID = 0; 699c7cc6ee6SMatthew G. Knepley PetscScalar *coordsIn = NULL; 700c7cc6ee6SMatthew G. Knepley PetscScalar *coords; 7013f6dc66eSMichael Lange PetscSection coordSection; 7023f6dc66eSMichael Lange Vec coordinates; 7035733454dSMatthew G. Knepley PetscInt coordSize, maxFaceSize = 0, totFaceVert = 0, f; 704c7cc6ee6SMatthew G. Knepley PetscMPIInt rank; 7052f0bd6dcSMichael Lange 7062f0bd6dcSMichael Lange PetscFunctionBegin; 7079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7082f0bd6dcSMichael Lange 709dd400576SPatrick Sanan if (rank == 0) { 7102f0bd6dcSMichael Lange FluentSection s; 7115733454dSMatthew G. Knepley 7125733454dSMatthew G. Knepley s.data = NULL; 713c7cc6ee6SMatthew G. Knepley numFaces = PETSC_DETERMINE; 7142f0bd6dcSMichael Lange do { 7159566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadSection(viewer, &s)); 7162f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 717f7320561SMichael Lange dim = s.nd; 718c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found dimension: %" PetscInt_FMT "\n", dim)); 71919d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 720c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 721c7cc6ee6SMatthew G. Knepley numVertices = s.last; 722c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of vertices: %" PetscInt_FMT "\n", numVertices)); 723c7cc6ee6SMatthew G. Knepley } else { 724c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found vertex coordinates\n")); 72528b400f6SJacob Faibussowitsch PetscCheck(!coordsIn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 7261a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *)s.data; 7273f6dc66eSMichael Lange } 7282f0bd6dcSMichael Lange 72919d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 730c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 731c7cc6ee6SMatthew G. Knepley numCells = s.last; 732c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cells %" PetscInt_FMT "\n", numCells)); 733c7cc6ee6SMatthew G. Knepley } else { 7345733454dSMatthew G. Knepley PetscCall(PetscMalloc2(numCells, &cellSizes, numCells, &cellTypes)); 735*9de2bd2cSJose E. Roman for (PetscInt c = 0; c < numCells; ++c) PetscCall(GetNumCellFaces(s.nd ? s.nd : (int)((PetscInt *)s.data)[c], &cellSizes[c], &cellTypes[c])); 7365733454dSMatthew G. Knepley PetscCall(PetscFree(s.data)); 7375733454dSMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cell faces %" PetscInt_FMT "\n", numCells && s.nd ? cellSizes[0] : 0)); 738f7320561SMichael Lange } 73919d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 740f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 741930bae4bSMatthew G. Knepley numFaces = (PetscInt)(s.last - s.first + 1); 7425733454dSMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of faces %" PetscInt_FMT " face vertices: %d\n", numFaces, s.nd)); 743f7320561SMichael Lange } else { /* Data section */ 7445733454dSMatthew G. Knepley PetscInt *tmp; 7455733454dSMatthew G. Knepley PetscInt totSize = 0, offset = 0, doffset; 746cf554e2cSMatthew G. Knepley 74708401ef6SPierre Jolivet PetscCheck(numFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 7485733454dSMatthew G. Knepley if (!faceZoneIDs) PetscCall(PetscMalloc3(numFaces, &faceSizes, numFaces * 2, &faceAdjCell, numFaces, &faceZoneIDs)); 7495733454dSMatthew G. Knepley // Record the zoneID and face size for each face set 7505733454dSMatthew G. Knepley for (unsigned int z = s.first - 1; z < s.last; z++) { 7515733454dSMatthew G. Knepley faceZoneIDs[z] = s.zoneID; 7525733454dSMatthew G. Knepley if (s.nd) { 7535733454dSMatthew G. Knepley faceSizes[z] = s.nd; 7545733454dSMatthew G. Knepley } else { 7555733454dSMatthew G. Knepley faceSizes[z] = ((PetscInt *)s.data)[offset]; 7565733454dSMatthew G. Knepley offset += faceSizes[z] + 3; 7575733454dSMatthew G. Knepley } 7585733454dSMatthew G. Knepley totSize += faceSizes[z]; 7595733454dSMatthew G. Knepley maxFaceSize = PetscMax(maxFaceSize, faceSizes[z]); 7605733454dSMatthew G. Knepley } 7615733454dSMatthew G. Knepley 7625733454dSMatthew G. Knepley offset = totFaceVert; 7635733454dSMatthew G. Knepley doffset = s.nd ? 0 : 1; 7645733454dSMatthew G. Knepley PetscCall(PetscMalloc1(totFaceVert + totSize, &tmp)); 7655733454dSMatthew G. Knepley if (faces) PetscCall(PetscArraycpy(tmp, faces, totFaceVert)); 7665733454dSMatthew G. Knepley PetscCall(PetscFree(faces)); 7675733454dSMatthew G. Knepley totFaceVert += totSize; 7685733454dSMatthew G. Knepley faces = tmp; 7695733454dSMatthew G. Knepley 7705733454dSMatthew G. Knepley // Record face vertices and adjacent faces 7715733454dSMatthew G. Knepley const PetscInt Nfz = s.last - s.first + 1; 7725733454dSMatthew G. Knepley for (PetscInt f = 0; f < Nfz; ++f) { 7735733454dSMatthew G. Knepley const PetscInt face = f + s.first - 1; 7745733454dSMatthew G. Knepley const PetscInt faceSize = faceSizes[face]; 7755733454dSMatthew G. Knepley 7765733454dSMatthew G. Knepley for (PetscInt v = 0; v < faceSize; ++v) faces[offset + v] = ((PetscInt *)s.data)[doffset + v]; 7775733454dSMatthew G. Knepley faceAdjCell[face * 2 + 0] = ((PetscInt *)s.data)[doffset + faceSize + 0]; 7785733454dSMatthew G. Knepley faceAdjCell[face * 2 + 1] = ((PetscInt *)s.data)[doffset + faceSize + 1]; 7795733454dSMatthew G. Knepley offset += faceSize; 7805733454dSMatthew G. Knepley doffset += faceSize + (s.nd ? 2 : 3); 7815733454dSMatthew G. Knepley } 7829566063dSJacob Faibussowitsch PetscCall(PetscFree(s.data)); 783f7320561SMichael Lange } 784c7cc6ee6SMatthew G. Knepley } else if (s.index == 39) { /* Label information */ 785c7cc6ee6SMatthew G. Knepley if (s.zoneID >= maxZoneID) { 786c7cc6ee6SMatthew G. Knepley DMLabel *tmpL; 787c7cc6ee6SMatthew G. Knepley const char **tmp; 788c7cc6ee6SMatthew G. Knepley unsigned int newmax = maxZoneID + 1; 789c7cc6ee6SMatthew G. Knepley 790c7cc6ee6SMatthew G. Knepley while (newmax < s.zoneID + 1) newmax *= 2; 791c7cc6ee6SMatthew G. Knepley PetscCall(PetscCalloc2(newmax, &tmp, newmax, &tmpL)); 792c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) { 793c7cc6ee6SMatthew G. Knepley tmp[i] = zoneNames[i]; 794c7cc6ee6SMatthew G. Knepley tmpL[i] = zoneLabels[i]; 795c7cc6ee6SMatthew G. Knepley } 796c7cc6ee6SMatthew G. Knepley maxZoneID = newmax; 797c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 798c7cc6ee6SMatthew G. Knepley zoneNames = tmp; 799c7cc6ee6SMatthew G. Knepley zoneLabels = tmpL; 800c7cc6ee6SMatthew G. Knepley } 801c7cc6ee6SMatthew G. Knepley zoneNames[s.zoneID] = (const char *)s.data; 8022f0bd6dcSMichael Lange } 8032f0bd6dcSMichael Lange } while (s.index >= 0); 8042f0bd6dcSMichael Lange } 8059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm)); 80608401ef6SPierre Jolivet PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 807f7320561SMichael Lange 808f7320561SMichael Lange /* Allocate cell-vertex mesh */ 8099566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 8109566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 8119566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 8125733454dSMatthew G. Knepley // We do not want this label automatically computed, instead we fill it here 8135733454dSMatthew G. Knepley PetscCall(DMCreateLabel(*dm, "celltype")); 814c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetChart(*dm, 0, numCells + numFaces + numVertices)); 815dd400576SPatrick Sanan if (rank == 0) { 81608401ef6SPierre Jolivet PetscCheck(numCells >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 8175733454dSMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) { 8185733454dSMatthew G. Knepley PetscCall(DMPlexSetConeSize(*dm, c, cellSizes[c])); 8195733454dSMatthew G. Knepley PetscCall(DMPlexSetCellType(*dm, c, cellTypes[c])); 820c7cc6ee6SMatthew G. Knepley } 8215733454dSMatthew G. Knepley for (PetscInt v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 8225733454dSMatthew G. Knepley for (PetscInt f = 0; f < numFaces; ++f) { 8235733454dSMatthew G. Knepley DMPolytopeType ct; 8245733454dSMatthew G. Knepley 8255733454dSMatthew G. Knepley switch (faceSizes[f]) { 8265733454dSMatthew G. Knepley case 2: 8275733454dSMatthew G. Knepley ct = DM_POLYTOPE_SEGMENT; 8285733454dSMatthew G. Knepley break; 8295733454dSMatthew G. Knepley case 3: 8305733454dSMatthew G. Knepley ct = DM_POLYTOPE_TRIANGLE; 8315733454dSMatthew G. Knepley break; 8325733454dSMatthew G. Knepley case 4: 8335733454dSMatthew G. Knepley ct = DM_POLYTOPE_QUADRILATERAL; 8345733454dSMatthew G. Knepley break; 8355733454dSMatthew G. Knepley default: 8365733454dSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file with cone size %" PetscInt_FMT, faceSizes[f]); 8375733454dSMatthew G. Knepley } 8385733454dSMatthew G. Knepley PetscCall(DMPlexSetConeSize(*dm, f + numCells + numVertices, faceSizes[f])); 8395733454dSMatthew G. Knepley PetscCall(DMPlexSetCellType(*dm, f + numCells + numVertices, ct)); 8405733454dSMatthew G. Knepley } 841f7320561SMichael Lange } 8429566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 843f7320561SMichael Lange 844dd400576SPatrick Sanan if (rank == 0 && faces) { 8455733454dSMatthew G. Knepley PetscSection s; 8465733454dSMatthew G. Knepley PetscInt *cones, csize, foffset = 0; 847f7320561SMichael Lange 848c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCones(*dm, &cones)); 8495733454dSMatthew G. Knepley PetscCall(DMPlexGetConeSection(*dm, &s)); 8505733454dSMatthew G. Knepley PetscCall(PetscSectionGetConstrainedStorageSize(s, &csize)); 8515733454dSMatthew G. Knepley for (PetscInt c = 0; c < csize; ++c) cones[c] = -1; 852c7cc6ee6SMatthew G. Knepley for (PetscInt f = 0; f < numFaces; f++) { 8535733454dSMatthew G. Knepley const PetscInt cl = faceAdjCell[f * 2 + 0] - 1; 8545733454dSMatthew G. Knepley const PetscInt cr = faceAdjCell[f * 2 + 1] - 1; 855c7cc6ee6SMatthew G. Knepley const PetscInt face = f + numCells + numVertices; 856c7cc6ee6SMatthew G. Knepley PetscInt fcone[16]; 857c7cc6ee6SMatthew G. Knepley 858c7cc6ee6SMatthew G. Knepley // How could Fluent define the outward normal differently? Is there no end to the pain? 859c7cc6ee6SMatthew G. Knepley if (dim == 3) { 860c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, -1)); 861c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, 0)); 862c7cc6ee6SMatthew G. Knepley } else { 863c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, 0)); 864c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, -1)); 8659371c9d4SSatish Balay } 8665733454dSMatthew G. Knepley PetscCheck(faceSizes[f] < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of face vertices %" PetscInt_FMT " exceeds temporary storage", faceSizes[f]); 8675733454dSMatthew G. Knepley for (PetscInt v = 0; v < faceSizes[f]; ++v) fcone[v] = faces[foffset + v] + numCells - 1; 8685733454dSMatthew G. Knepley foffset += faceSizes[f]; 869c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(*dm, face, fcone)); 870f7320561SMichael Lange } 871f7320561SMichael Lange } 8729566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 8739566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 874c7cc6ee6SMatthew G. Knepley if (dim == 3) { 8755fd9971aSMatthew G. Knepley DM idm; 876f7320561SMichael Lange 877c7cc6ee6SMatthew G. Knepley PetscCall(DMCreate(PetscObjectComm((PetscObject)*dm), &idm)); 878c7cc6ee6SMatthew G. Knepley PetscCall(DMSetType(idm, DMPLEX)); 879c7cc6ee6SMatthew G. Knepley PetscCall(DMSetDimension(idm, dim)); 880c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInterpolateFaces_Internal(*dm, 1, idm)); 8819566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 882f7320561SMichael Lange *dm = idm; 883f7320561SMichael Lange } 884c7cc6ee6SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-cas_dm_view")); 885c7cc6ee6SMatthew G. Knepley if (rank == 0 && faces) { 8865733454dSMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) PetscCall(ReorderCell(viewer, *dm, c, cellTypes[c])); 887c7cc6ee6SMatthew G. Knepley } 888f7320561SMichael Lange 889dd400576SPatrick Sanan if (rank == 0 && faces) { 8905733454dSMatthew G. Knepley PetscInt joinSize, meetSize, *fverts, cells[2]; 891631eb916SMichael Lange const PetscInt *join, *meet; 8925733454dSMatthew G. Knepley PetscInt foffset = 0; 8935733454dSMatthew G. Knepley 8945733454dSMatthew G. Knepley PetscCall(PetscMalloc1(maxFaceSize, &fverts)); 895ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 896ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 8975733454dSMatthew G. Knepley const PetscInt cl = faceAdjCell[f * 2 + 0] - 1; 8985733454dSMatthew G. Knepley const PetscInt cr = faceAdjCell[f * 2 + 1] - 1; 899c7cc6ee6SMatthew G. Knepley const PetscInt id = (PetscInt)faceZoneIDs[f]; 900c7cc6ee6SMatthew G. Knepley 901631eb916SMichael Lange if (cl > 0 && cr > 0) { 902631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 9039371c9d4SSatish Balay cells[0] = cl; 9049371c9d4SSatish Balay cells[1] = cr; 9059566063dSJacob Faibussowitsch PetscCall(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet)); 906c7cc6ee6SMatthew G. Knepley PetscCheck(meetSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT " cells: %" PetscInt_FMT ", %" PetscInt_FMT, f, cl, cr); 907c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], id)); 908c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], meet[0], 1)); 9095733454dSMatthew G. Knepley PetscCall(DMPlexRestoreMeet(*dm, meetSize, fverts, &meetSize, &meet)); 910631eb916SMichael Lange } else { 9115733454dSMatthew G. Knepley for (PetscInt fi = 0; fi < faceSizes[f]; fi++) fverts[fi] = faces[foffset + fi] + numCells - 1; 9125733454dSMatthew G. Knepley PetscCall(DMPlexGetFullJoin(*dm, faceSizes[f], fverts, &joinSize, &join)); 91363a3b9bcSJacob Faibussowitsch PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f); 914c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], id)); 915c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], join[0], 1)); 9165733454dSMatthew G. Knepley PetscCall(DMPlexRestoreJoin(*dm, joinSize, fverts, &joinSize, &join)); 917ec78a56aSMichael Lange } 9185733454dSMatthew G. Knepley foffset += faceSizes[f]; 919631eb916SMichael Lange } 9209566063dSJacob Faibussowitsch PetscCall(PetscFree(fverts)); 921ec78a56aSMichael Lange } 922ec78a56aSMichael Lange 9237368db69SLisandro Dalcin { /* Create Face Sets label at all processes */ 9249371c9d4SSatish Balay enum { 9259371c9d4SSatish Balay n = 1 9269371c9d4SSatish Balay }; 9277368db69SLisandro Dalcin PetscBool flag[n]; 9287368db69SLisandro Dalcin 9297368db69SLisandro Dalcin flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE; 9309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 9319566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 932c7cc6ee6SMatthew G. Knepley // TODO Code to create all the zone labels on each process 933c7cc6ee6SMatthew G. Knepley } 934c7cc6ee6SMatthew G. Knepley 935c7cc6ee6SMatthew G. Knepley if (!interpolate) { 936c7cc6ee6SMatthew G. Knepley DM udm; 937c7cc6ee6SMatthew G. Knepley 938c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexUninterpolate(*dm, &udm)); 939c7cc6ee6SMatthew G. Knepley PetscCall(DMDestroy(dm)); 940c7cc6ee6SMatthew G. Knepley *dm = udm; 9417368db69SLisandro Dalcin } 9427368db69SLisandro Dalcin 9433f6dc66eSMichael Lange /* Read coordinates */ 9449566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 9459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 9469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 9479566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 948c7cc6ee6SMatthew G. Knepley for (PetscInt v = numCells; v < numCells + numVertices; ++v) { 9499566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dim)); 9509566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 9513f6dc66eSMichael Lange } 9529566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 9539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 9549566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 9559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 9569566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 9579566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 9589566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 959dd400576SPatrick Sanan if (rank == 0 && coordsIn) { 960c7cc6ee6SMatthew G. Knepley for (PetscInt v = 0; v < numVertices; ++v) { 961c7cc6ee6SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) coords[v * dim + d] = coordsIn[v * dim + d]; 9623f6dc66eSMichael Lange } 9633f6dc66eSMichael Lange } 9649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 9659566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 9669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 9677368db69SLisandro Dalcin 968dd400576SPatrick Sanan if (rank == 0) { 9699566063dSJacob Faibussowitsch PetscCall(PetscFree(cellVertices)); 9705733454dSMatthew G. Knepley PetscCall(PetscFree2(cellSizes, cellTypes)); 9719566063dSJacob Faibussowitsch PetscCall(PetscFree(faces)); 9725733454dSMatthew G. Knepley PetscCall(PetscFree3(faceSizes, faceAdjCell, faceZoneIDs)); 9739566063dSJacob Faibussowitsch PetscCall(PetscFree(coordsIn)); 974c7cc6ee6SMatthew G. Knepley if (zoneNames) 975c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) PetscCall(PetscFree(zoneNames[i])); 976c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 977f7320561SMichael Lange } 9783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9792f0bd6dcSMichael Lange } 980