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 5cc4c1da9SBarry Smith /*@ 6a1cb98faSBarry Smith DMPlexCreateFluentFromFile - Create a `DMPLEX` mesh from a Fluent mesh file 72f0bd6dcSMichael Lange 8a1cb98faSBarry Smith Collective 9a1cb98faSBarry Smith 10a1cb98faSBarry Smith Input Parameters: 112f0bd6dcSMichael Lange + comm - The MPI communicator 122f0bd6dcSMichael Lange . filename - Name of the Fluent mesh file 132f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 142f0bd6dcSMichael Lange 152f0bd6dcSMichael Lange Output Parameter: 16a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 172f0bd6dcSMichael Lange 182f0bd6dcSMichael Lange Level: beginner 192f0bd6dcSMichael Lange 201cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMPlexCreateFromFile()`, `DMPlexCreateFluent()`, `DMPlexCreate()` 212f0bd6dcSMichael Lange @*/ 22d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 23d71ae5a4SJacob Faibussowitsch { 242f0bd6dcSMichael Lange PetscViewer viewer; 252f0bd6dcSMichael Lange 262f0bd6dcSMichael Lange PetscFunctionBegin; 272f0bd6dcSMichael Lange /* Create file viewer and build plex */ 289566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, &viewer)); 299566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII)); 309566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ)); 319566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(viewer, filename)); 329566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent(comm, viewer, interpolate, dm)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 352f0bd6dcSMichael Lange } 362f0bd6dcSMichael Lange 37d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 38d71ae5a4SJacob Faibussowitsch { 394b375a39SMichael Lange PetscInt ret, i = 0; 402f0bd6dcSMichael Lange 412f0bd6dcSMichael Lange PetscFunctionBegin; 42f4f49eeaSPierre Jolivet do PetscCall(PetscViewerRead(viewer, &buffer[i++], 1, &ret, PETSC_CHAR)); 43c7cc6ee6SMatthew G. Knepley while (ret > 0 && buffer[i - 1] != '\0' && buffer[i - 1] != delim && i < PETSC_MAX_PATH_LEN - 1); 449371c9d4SSatish Balay if (!ret) buffer[i - 1] = '\0'; 459371c9d4SSatish Balay else buffer[i] = '\0'; 46c7cc6ee6SMatthew 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."); 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 482f0bd6dcSMichael Lange } 492f0bd6dcSMichael Lange 50c7cc6ee6SMatthew G. Knepley static PetscErrorCode DMPlexCreateFluent_ReadValues(PetscViewer viewer, void *data, PetscInt count, PetscDataType dtype, PetscBool binary, PetscInt *numClosingParens) 51d71ae5a4SJacob Faibussowitsch { 5270c9a859SSatish Balay int fdes = 0; 5319d58f9dSMichael Lange FILE *file; 5419d58f9dSMichael Lange PetscInt i; 5519d58f9dSMichael Lange 5619d58f9dSMichael Lange PetscFunctionBegin; 57c7cc6ee6SMatthew G. Knepley *numClosingParens = 0; 5819d58f9dSMichael Lange if (binary) { 5919d58f9dSMichael Lange /* Extract raw file descriptor to read binary block */ 609566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIGetPointer(viewer, &file)); 61c69effb2SJacob Faibussowitsch PetscCall(PetscFFlush(file)); 629371c9d4SSatish Balay fdes = fileno(file); 6319d58f9dSMichael Lange } 6419d58f9dSMichael Lange 6519d58f9dSMichael Lange if (!binary && dtype == PETSC_INT) { 6619d58f9dSMichael Lange char cbuf[256]; 67cfb60857SMatthew G. Knepley unsigned int ibuf; 68cfb60857SMatthew G. Knepley int snum; 6919d58f9dSMichael Lange /* Parse hexadecimal ascii integers */ 7019d58f9dSMichael Lange for (i = 0; i < count; i++) { 71c7cc6ee6SMatthew G. Knepley size_t len; 72c7cc6ee6SMatthew G. Knepley 739566063dSJacob Faibussowitsch PetscCall(PetscViewerRead(viewer, cbuf, 1, NULL, PETSC_STRING)); 7419d58f9dSMichael Lange snum = sscanf(cbuf, "%x", &ibuf); 7508401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 7619d58f9dSMichael Lange ((PetscInt *)data)[i] = (PetscInt)ibuf; 77c7cc6ee6SMatthew G. Knepley // Check for trailing parentheses 78c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrlen(cbuf, &len)); 79c7cc6ee6SMatthew G. Knepley while (cbuf[len - 1] == ')' && len > 0) { 80c7cc6ee6SMatthew G. Knepley ++(*numClosingParens); 81c7cc6ee6SMatthew G. Knepley --len; 82c7cc6ee6SMatthew G. Knepley } 8319d58f9dSMichael Lange } 8419d58f9dSMichael Lange } else if (binary && dtype == PETSC_INT) { 8519d58f9dSMichael Lange /* Always read 32-bit ints and cast to PetscInt */ 8619d58f9dSMichael Lange int *ibuf; 879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, &ibuf)); 889566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes, ibuf, count, NULL, PETSC_ENUM)); 899566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(ibuf, PETSC_ENUM, count)); 90835f2295SStefano Zampini for (i = 0; i < count; i++) ((PetscInt *)data)[i] = ibuf[i]; 919566063dSJacob Faibussowitsch PetscCall(PetscFree(ibuf)); 9219d58f9dSMichael Lange 9319d58f9dSMichael Lange } else if (binary && dtype == PETSC_SCALAR) { 9419d58f9dSMichael Lange float *fbuf; 9519d58f9dSMichael Lange /* Always read 32-bit floats and cast to PetscScalar */ 969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(count, &fbuf)); 979566063dSJacob Faibussowitsch PetscCall(PetscBinaryRead(fdes, fbuf, count, NULL, PETSC_FLOAT)); 989566063dSJacob Faibussowitsch PetscCall(PetscByteSwap(fbuf, PETSC_FLOAT, count)); 99835f2295SStefano Zampini for (i = 0; i < count; i++) ((PetscScalar *)data)[i] = fbuf[i]; 1009566063dSJacob Faibussowitsch PetscCall(PetscFree(fbuf)); 10119d58f9dSMichael Lange } else { 1029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIRead(viewer, data, count, NULL, dtype)); 10319d58f9dSMichael Lange } 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10519d58f9dSMichael Lange } 10619d58f9dSMichael Lange 107d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 108d71ae5a4SJacob Faibussowitsch { 1092f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 1102f0bd6dcSMichael Lange int snum; 1112f0bd6dcSMichael Lange 1122f0bd6dcSMichael Lange PetscFunctionBegin; 1132f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 1149566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1159566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ' ')); 116f4f49eeaSPierre Jolivet snum = sscanf(buffer, "%d", &s->index); 1172f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 1189371c9d4SSatish Balay if (snum < 1) { 1199371c9d4SSatish Balay s->index = -1; 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1219371c9d4SSatish Balay } 1222f0bd6dcSMichael Lange 1232f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 1249566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 1252f0bd6dcSMichael Lange 1262f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 1279566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 128f4f49eeaSPierre Jolivet snum = sscanf(buffer, "%d", &s->nd); 12908401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1302f0bd6dcSMichael Lange 13119d58f9dSMichael Lange } else if (s->index == 10 || s->index == 2010) { /* Vertices */ 132c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 133c7cc6ee6SMatthew G. Knepley 1349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 135f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 13608401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1373f6dc66eSMichael Lange if (s->zoneID > 0) { 1383f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 1399566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(s->nd * numCoords, (PetscScalar **)&s->data)); 141c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd * numCoords, PETSC_SCALAR, s->index == 2010 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 142c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 143c7cc6ee6SMatthew G. Knepley else --numClosingParens; 1443f6dc66eSMichael Lange } 145c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 146c7cc6ee6SMatthew G. Knepley else --numClosingParens; 147c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 14819d58f9dSMichael Lange } else if (s->index == 12 || s->index == 2012) { /* Cells */ 149c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 150c7cc6ee6SMatthew G. Knepley 1519566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 152f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x", &s->zoneID); 15308401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1542f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 155f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd); 15608401ef6SPierre Jolivet PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 157f7320561SMichael Lange } else { /* Data section */ 158f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 15908401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 160895fcc3eSMichael Lange if (s->nd == 0) { 161895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 16219d58f9dSMichael Lange PetscInt numCells = s->last - s->first + 1; 1639566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 1649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCells, (PetscInt **)&s->data)); 165c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, numCells, PETSC_INT, s->index == 2012 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(s->data)); 167c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 168c7cc6ee6SMatthew G. Knepley else --numClosingParens; 169895fcc3eSMichael Lange } 1702f0bd6dcSMichael Lange } 171c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 172c7cc6ee6SMatthew G. Knepley else --numClosingParens; 173c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 17419d58f9dSMichael Lange } else if (s->index == 13 || s->index == 2013) { /* Faces */ 175c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 176c7cc6ee6SMatthew G. Knepley 1779566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 178f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x", &s->zoneID); 17908401ef6SPierre Jolivet PetscCheck(snum == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1802f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 181f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d)", &s->zoneID, &s->first, &s->last, &s->nd); 18208401ef6SPierre Jolivet PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 183f7320561SMichael Lange } else { /* Data section */ 184137d0469SJed Brown PetscInt f, numEntries, numFaces; 185f4f49eeaSPierre Jolivet snum = sscanf(buffer, "(%x %x %x %d %d)", &s->zoneID, &s->first, &s->last, &s->type, &s->nd); 18608401ef6SPierre Jolivet PetscCheck(snum == 5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1879566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '(')); 188f7320561SMichael Lange switch (s->nd) { 189d71ae5a4SJacob Faibussowitsch case 0: 190d71ae5a4SJacob Faibussowitsch numEntries = PETSC_DETERMINE; 191d71ae5a4SJacob Faibussowitsch break; 192d71ae5a4SJacob Faibussowitsch case 2: 193d71ae5a4SJacob Faibussowitsch numEntries = 2 + 2; 194d71ae5a4SJacob Faibussowitsch break; /* linear */ 195d71ae5a4SJacob Faibussowitsch case 3: 196d71ae5a4SJacob Faibussowitsch numEntries = 2 + 3; 197d71ae5a4SJacob Faibussowitsch break; /* triangular */ 198d71ae5a4SJacob Faibussowitsch case 4: 199d71ae5a4SJacob Faibussowitsch numEntries = 2 + 4; 200d71ae5a4SJacob Faibussowitsch break; /* quadrilateral */ 201d71ae5a4SJacob Faibussowitsch default: 202d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 203f7320561SMichael Lange } 204f7320561SMichael Lange numFaces = s->last - s->first + 1; 205895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 206895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 2079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data)); 208895fcc3eSMichael Lange } 209f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 210895fcc3eSMichael Lange if (s->nd == 0) { 211895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 212137d0469SJed Brown PetscInt numFaceVert = 0; 213c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, &numFaceVert, 1, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 214895fcc3eSMichael Lange if (numEntries == PETSC_DETERMINE) { 215895fcc3eSMichael Lange numEntries = numFaceVert + 2; 2169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numEntries * numFaces, (PetscInt **)&s->data)); 217895fcc3eSMichael Lange } else { 2181dca8a05SBarry Smith PetscCheck(numEntries == numFaceVert + 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 219895fcc3eSMichael Lange } 220895fcc3eSMichael Lange } 221c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, &(((PetscInt *)s->data)[f * numEntries]), numEntries, PETSC_INT, s->index == 2013 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 222f7320561SMichael Lange } 223835f2295SStefano Zampini PetscCall(PetscMPIIntCast(numEntries - 2, &s->nd)); 224c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 225c7cc6ee6SMatthew G. Knepley else --numClosingParens; 2262f0bd6dcSMichael Lange } 227c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 228c7cc6ee6SMatthew G. Knepley else --numClosingParens; 229c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 230c7cc6ee6SMatthew G. Knepley } else if (s->index == 39) { /* Label information */ 231c7cc6ee6SMatthew G. Knepley char labelName[PETSC_MAX_PATH_LEN]; 232c7cc6ee6SMatthew G. Knepley char caseName[PETSC_MAX_PATH_LEN]; 2332f0bd6dcSMichael Lange 234c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 235c3279546SPierre Jolivet snum = sscanf(buffer, "(%u %s %s %d)", &s->zoneID, caseName, labelName, &s->nd); 236c7cc6ee6SMatthew G. Knepley PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file: %d", snum); 237c7cc6ee6SMatthew G. Knepley PetscInt depth = 1; 238c7cc6ee6SMatthew G. Knepley do { 239c7cc6ee6SMatthew G. Knepley /* Match parentheses when parsing unknown sections */ 240c7cc6ee6SMatthew G. Knepley do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR)); 241c7cc6ee6SMatthew G. Knepley while (buffer[0] != '(' && buffer[0] != ')'); 242c7cc6ee6SMatthew G. Knepley if (buffer[0] == '(') depth++; 243c7cc6ee6SMatthew G. Knepley if (buffer[0] == ')') depth--; 244c7cc6ee6SMatthew G. Knepley } while (depth > 0); 245c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 246c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrallocpy(labelName, (char **)&s->data)); 247c3279546SPierre Jolivet PetscCall(PetscInfo((PetscObject)viewer, "CASE: Zone ID %u is label %s\n", s->zoneID, labelName)); 2482f0bd6dcSMichael Lange } else { /* Unknown section type */ 249060da220SMatthew G. Knepley PetscInt depth = 1; 2502f0bd6dcSMichael Lange do { 2512f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 252f4f49eeaSPierre Jolivet do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR)); 2532f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 2542f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 2552f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 2562f0bd6dcSMichael Lange } while (depth > 0); 2579566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 2582f0bd6dcSMichael Lange } 2593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2602f0bd6dcSMichael Lange } 2612f0bd6dcSMichael Lange 262c7cc6ee6SMatthew G. Knepley // Inserts point `face` with orientation `ornt` into the cone of point `cell` at position `c`, which is the first empty slot 263c7cc6ee6SMatthew G. Knepley static PetscErrorCode InsertFace(DM dm, PetscInt cell, PetscInt face, PetscInt ornt) 264c7cc6ee6SMatthew G. Knepley { 265c7cc6ee6SMatthew G. Knepley const PetscInt *cone; 266c7cc6ee6SMatthew G. Knepley PetscInt coneSize, c; 267c7cc6ee6SMatthew G. Knepley 268c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 269c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cell, &cone)); 270c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 271c7cc6ee6SMatthew G. Knepley for (c = 0; c < coneSize; ++c) 272c7cc6ee6SMatthew G. Knepley if (cone[c] < 0) break; 273c7cc6ee6SMatthew 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, face, cell); 274c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertCone(dm, cell, c, face)); 275c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertConeOrientation(dm, cell, c, ornt)); 276c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 277c7cc6ee6SMatthew G. Knepley } 278c7cc6ee6SMatthew G. Knepley 279c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderPolygon(DM dm, PetscInt cell) 280c7cc6ee6SMatthew G. Knepley { 281c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt; 282c7cc6ee6SMatthew G. Knepley PetscInt coneSize, newCone[16], newOrnt[16]; 283c7cc6ee6SMatthew G. Knepley 284c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 285c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 286c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 287c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 288c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 289c7cc6ee6SMatthew G. Knepley for (PetscInt c = 1; c < coneSize; ++c) { 290c7cc6ee6SMatthew G. Knepley const PetscInt *fcone; 291c7cc6ee6SMatthew G. Knepley PetscInt firstVertex, lastVertex, c2; 292c7cc6ee6SMatthew G. Knepley 293c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[c - 1], &fcone)); 294c7cc6ee6SMatthew G. Knepley lastVertex = newOrnt[c - 1] ? fcone[0] : fcone[1]; 295c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < coneSize; ++c2) { 296c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2; 297c7cc6ee6SMatthew G. Knepley 298c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cone[c2], &fcone2)); 299c7cc6ee6SMatthew G. Knepley firstVertex = ornt[c2] ? fcone2[1] : fcone2[0]; 300c7cc6ee6SMatthew G. Knepley if (lastVertex == firstVertex) { 301c7cc6ee6SMatthew G. Knepley // Point `cell` matched point `lastVertex` on face `cone[c2]` with orientation `ornt[c2]` 302c7cc6ee6SMatthew G. Knepley break; 303c7cc6ee6SMatthew G. Knepley } 304c7cc6ee6SMatthew G. Knepley } 305c7cc6ee6SMatthew 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); 306c7cc6ee6SMatthew G. Knepley newCone[c] = cone[c2]; 307c7cc6ee6SMatthew G. Knepley newOrnt[c] = ornt[c2]; 308c7cc6ee6SMatthew G. Knepley } 309c7cc6ee6SMatthew G. Knepley { 310c7cc6ee6SMatthew G. Knepley const PetscInt *fcone, *fcone2; 311c7cc6ee6SMatthew G. Knepley PetscInt vertex, vertex2; 312c7cc6ee6SMatthew G. Knepley 313c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[coneSize - 1], &fcone)); 314c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[0], &fcone2)); 315c7cc6ee6SMatthew G. Knepley vertex = newOrnt[coneSize - 1] ? fcone[0] : fcone[1]; 316c7cc6ee6SMatthew G. Knepley vertex2 = newOrnt[0] ? fcone2[1] : fcone2[0]; 317c7cc6ee6SMatthew G. Knepley PetscCheck(vertex == vertex2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " did not match at the endpoint", cell); 318c7cc6ee6SMatthew G. Knepley } 319c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 320c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 321c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 322c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 323c7cc6ee6SMatthew G. Knepley } 324c7cc6ee6SMatthew G. Knepley 325*fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderTetrahedron(PetscViewer viewer, DM dm, PetscInt cell) 326c7cc6ee6SMatthew G. Knepley { 327c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr, faces[4] = {0, 1, 3, 2}; 328c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 329c7cc6ee6SMatthew G. Knepley 330c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 331c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 332c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 333c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 334c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 335c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]); 336c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial triangle 337c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 3; ++e) { 338c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 339c7cc6ee6SMatthew G. Knepley PetscInt c; 340c7cc6ee6SMatthew G. Knepley 341c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the tetrahedron 342c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 343c7cc6ee6SMatthew G. Knepley for (c = 1; c < 4; ++c) { 344c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 345c7cc6ee6SMatthew G. Knepley PetscInt c2; 346*fc8f8d65SMatthew G. Knepley PetscBool flip = PETSC_FALSE; 347c7cc6ee6SMatthew G. Knepley 348c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 349c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 350c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]); 351c7cc6ee6SMatthew G. Knepley // Check for edge 352c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 3; ++c2) { 353*fc8f8d65SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 354c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 355c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 356*fc8f8d65SMatthew 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); 357c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 358*fc8f8d65SMatthew 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)); 359*fc8f8d65SMatthew G. Knepley flip = eornt != -(eornt2 + 1) ? PETSC_TRUE : PETSC_FALSE; 360c7cc6ee6SMatthew G. Knepley break; 361c7cc6ee6SMatthew G. Knepley } 362c7cc6ee6SMatthew G. Knepley } 363c7cc6ee6SMatthew G. Knepley if (c2 < 3) { 364c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 365c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 366*fc8f8d65SMatthew G. Knepley // Face 1 should match its edge 2 367*fc8f8d65SMatthew G. Knepley // Face 2 should match its edge 0 368*fc8f8d65SMatthew G. Knepley // Face 3 should match its edge 0 369*fc8f8d65SMatthew G. Knepley if (flip) { 370*fc8f8d65SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, -((c2 + (!e ? 1 : 2)) % 3 + 1), ornt[c]); 371*fc8f8d65SMatthew G. Knepley } else { 372c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, !e ? (c2 + 1) % 3 : c2, ornt[c]); 373*fc8f8d65SMatthew G. Knepley } 374c7cc6ee6SMatthew G. Knepley break; 375c7cc6ee6SMatthew G. Knepley } 376c7cc6ee6SMatthew G. Knepley } 377c7cc6ee6SMatthew 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); 378c7cc6ee6SMatthew G. Knepley } 379c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 380c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 381c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 382c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 383c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 384c7cc6ee6SMatthew G. Knepley } 385c7cc6ee6SMatthew G. Knepley 386c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderHexahedron(DM dm, PetscInt cell) 387c7cc6ee6SMatthew G. Knepley { 388c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr; 389c7cc6ee6SMatthew G. Knepley const PetscInt faces[6] = {0, 5, 3, 4, 2, 1}; 390c7cc6ee6SMatthew G. Knepley PetscInt used[6] = {1, 0, 0, 0, 0, 0}; 391c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 392c7cc6ee6SMatthew G. Knepley 393c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 394c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 395c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 396c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 397c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 398c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[0]); 399c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial quadrilateral 400c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 4; ++e) { 401c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 402c7cc6ee6SMatthew G. Knepley PetscInt c; 403c7cc6ee6SMatthew G. Knepley 404c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the hexahedron 405c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 406c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) { 407c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 408c7cc6ee6SMatthew G. Knepley PetscInt c2; 409c7cc6ee6SMatthew G. Knepley 410c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 411c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 412c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 413c7cc6ee6SMatthew G. Knepley // Check for edge 414c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 415c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 416c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 417c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 418c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 419c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 420c7cc6ee6SMatthew G. Knepley break; 421c7cc6ee6SMatthew G. Knepley } 422c7cc6ee6SMatthew G. Knepley } 423c7cc6ee6SMatthew G. Knepley if (c2 < 4) { 424c7cc6ee6SMatthew G. Knepley used[c] = 1; 425c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 426c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 427c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, !e ? (c2 + 1) % 4 : c2, ornt[c]); 428c7cc6ee6SMatthew G. Knepley break; 429c7cc6ee6SMatthew G. Knepley } 430c7cc6ee6SMatthew G. Knepley } 431c7cc6ee6SMatthew 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); 432c7cc6ee6SMatthew G. Knepley } 433c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 434c7cc6ee6SMatthew G. Knepley // Add last face 435c7cc6ee6SMatthew G. Knepley { 436c7cc6ee6SMatthew G. Knepley PetscInt c, c2; 437c7cc6ee6SMatthew G. Knepley 438c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) 439c7cc6ee6SMatthew G. Knepley if (!used[c]) break; 440c7cc6ee6SMatthew G. Knepley PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell); 441c7cc6ee6SMatthew G. Knepley // Match first edge to 3rd edge in newCone[2] 442c7cc6ee6SMatthew G. Knepley { 443c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 444c7cc6ee6SMatthew G. Knepley 445c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt)); 446c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]); 447c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 448c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 449c7cc6ee6SMatthew G. Knepley 450c7cc6ee6SMatthew G. Knepley const PetscInt e = 2; 451c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 452c7cc6ee6SMatthew 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]` 453c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 454c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 455c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 456c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 457c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 458c7cc6ee6SMatthew G. Knepley // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 459c7cc6ee6SMatthew G. Knepley break; 460c7cc6ee6SMatthew G. Knepley } 461c7cc6ee6SMatthew G. Knepley } 462c7cc6ee6SMatthew G. Knepley PetscCheck(c2 < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in"); 463c7cc6ee6SMatthew G. Knepley } 464c7cc6ee6SMatthew G. Knepley newCone[faces[5]] = cone[c]; 465c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched 466c7cc6ee6SMatthew G. Knepley newOrnt[faces[5]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]); 467c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 468c7cc6ee6SMatthew G. Knepley } 469c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 470c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 471c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 472c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 473c7cc6ee6SMatthew G. Knepley } 474c7cc6ee6SMatthew G. Knepley 475*fc8f8d65SMatthew G. Knepley static PetscErrorCode ReorderCell(PetscViewer viewer, DM dm, PetscInt cell, DMPolytopeType ct) 476c7cc6ee6SMatthew G. Knepley { 477c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 478c7cc6ee6SMatthew G. Knepley switch (ct) { 479c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 480c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 481c7cc6ee6SMatthew G. Knepley PetscCall(ReorderPolygon(dm, cell)); 482c7cc6ee6SMatthew G. Knepley break; 483c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TETRAHEDRON: 484*fc8f8d65SMatthew G. Knepley PetscCall(ReorderTetrahedron(viewer, dm, cell)); 485c7cc6ee6SMatthew G. Knepley break; 486c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_HEXAHEDRON: 487c7cc6ee6SMatthew G. Knepley PetscCall(ReorderHexahedron(dm, cell)); 488c7cc6ee6SMatthew G. Knepley break; 489c7cc6ee6SMatthew G. Knepley default: 490c7cc6ee6SMatthew G. Knepley PetscCheck(0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Celltype %s is unsupported", DMPolytopeTypes[ct]); 491c7cc6ee6SMatthew G. Knepley break; 492c7cc6ee6SMatthew G. Knepley } 493c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 494c7cc6ee6SMatthew G. Knepley } 495c7cc6ee6SMatthew G. Knepley 4962f0bd6dcSMichael Lange /*@C 4971d27aa22SBarry Smith DMPlexCreateFluent - Create a `DMPLEX` mesh from a Fluent mesh file <http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm>. 4982f0bd6dcSMichael Lange 499d083f849SBarry Smith Collective 5002f0bd6dcSMichael Lange 5012f0bd6dcSMichael Lange Input Parameters: 5022f0bd6dcSMichael Lange + comm - The MPI communicator 503a1cb98faSBarry Smith . viewer - The `PetscViewer` associated with a Fluent mesh file 5042f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 5052f0bd6dcSMichael Lange 5062f0bd6dcSMichael Lange Output Parameter: 507a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 5082f0bd6dcSMichael Lange 5092f0bd6dcSMichael Lange Level: beginner 5102f0bd6dcSMichael Lange 5111cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()` 5122f0bd6dcSMichael Lange @*/ 513d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 514d71ae5a4SJacob Faibussowitsch { 515c7cc6ee6SMatthew G. Knepley PetscInt dim = PETSC_DETERMINE; 516c7cc6ee6SMatthew G. Knepley PetscInt numCells = 0; 517c7cc6ee6SMatthew G. Knepley PetscInt numVertices = 0; 518c7cc6ee6SMatthew G. Knepley PetscInt numCellFaces = PETSC_DETERMINE; 519c7cc6ee6SMatthew G. Knepley DMPolytopeType ct = DM_POLYTOPE_UNKNOWN_CELL; 520c7cc6ee6SMatthew G. Knepley PetscInt numFaces = 0; 521c7cc6ee6SMatthew G. Knepley PetscInt numFaceEntries = PETSC_DETERMINE; 522c7cc6ee6SMatthew G. Knepley PetscInt numFaceVertices = PETSC_DETERMINE; 523c7cc6ee6SMatthew G. Knepley PetscInt *faces = NULL; 524c7cc6ee6SMatthew G. Knepley PetscInt *cellVertices = NULL; 525c7cc6ee6SMatthew G. Knepley unsigned int *faceZoneIDs = NULL; 5267368db69SLisandro Dalcin DMLabel faceSets = NULL; 527c7cc6ee6SMatthew G. Knepley DMLabel *zoneLabels = NULL; 528c7cc6ee6SMatthew G. Knepley const char **zoneNames = NULL; 529c7cc6ee6SMatthew G. Knepley unsigned int maxZoneID = 0; 530c7cc6ee6SMatthew G. Knepley PetscScalar *coordsIn = NULL; 531c7cc6ee6SMatthew G. Knepley PetscScalar *coords; 5323f6dc66eSMichael Lange PetscSection coordSection; 5333f6dc66eSMichael Lange Vec coordinates; 534c7cc6ee6SMatthew G. Knepley PetscInt coordSize, f; 535c7cc6ee6SMatthew G. Knepley PetscMPIInt rank; 5362f0bd6dcSMichael Lange 5372f0bd6dcSMichael Lange PetscFunctionBegin; 5389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5392f0bd6dcSMichael Lange 540dd400576SPatrick Sanan if (rank == 0) { 5412f0bd6dcSMichael Lange FluentSection s; 542c7cc6ee6SMatthew G. Knepley numFaces = PETSC_DETERMINE; 5432f0bd6dcSMichael Lange do { 5449566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadSection(viewer, &s)); 5452f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 546f7320561SMichael Lange dim = s.nd; 547c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found dimension: %" PetscInt_FMT "\n", dim)); 54819d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 549c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 550c7cc6ee6SMatthew G. Knepley numVertices = s.last; 551c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of vertices: %" PetscInt_FMT "\n", numVertices)); 552c7cc6ee6SMatthew G. Knepley } else { 553c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found vertex coordinates\n")); 55428b400f6SJacob Faibussowitsch PetscCheck(!coordsIn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 5551a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *)s.data; 5563f6dc66eSMichael Lange } 5572f0bd6dcSMichael Lange 55819d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 559c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 560c7cc6ee6SMatthew G. Knepley numCells = s.last; 561c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cells %" PetscInt_FMT "\n", numCells)); 562c7cc6ee6SMatthew G. Knepley } else { 563f7320561SMichael Lange switch (s.nd) { 564d71ae5a4SJacob Faibussowitsch case 0: 565c7cc6ee6SMatthew G. Knepley numCellFaces = PETSC_DETERMINE; 566c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_POINT; 567d71ae5a4SJacob Faibussowitsch break; 568d71ae5a4SJacob Faibussowitsch case 1: 569c7cc6ee6SMatthew G. Knepley numCellFaces = 3; 570c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TRIANGLE; 571c7cc6ee6SMatthew G. Knepley break; 572d71ae5a4SJacob Faibussowitsch case 2: 573c7cc6ee6SMatthew G. Knepley numCellFaces = 4; 574c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TETRAHEDRON; 575c7cc6ee6SMatthew G. Knepley break; 576d71ae5a4SJacob Faibussowitsch case 3: 577c7cc6ee6SMatthew G. Knepley numCellFaces = 4; 578c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_QUADRILATERAL; 579c7cc6ee6SMatthew G. Knepley break; 580d71ae5a4SJacob Faibussowitsch case 4: 581c7cc6ee6SMatthew G. Knepley numCellFaces = 6; 582c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_HEXAHEDRON; 583c7cc6ee6SMatthew G. Knepley break; 584d71ae5a4SJacob Faibussowitsch case 5: 585c7cc6ee6SMatthew G. Knepley numCellFaces = 5; 586c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_PYRAMID; 587c7cc6ee6SMatthew G. Knepley break; 588d71ae5a4SJacob Faibussowitsch case 6: 589c7cc6ee6SMatthew G. Knepley numCellFaces = 5; 590c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TRI_PRISM; 591c7cc6ee6SMatthew G. Knepley break; 592d71ae5a4SJacob Faibussowitsch default: 593c7cc6ee6SMatthew G. Knepley numCellFaces = PETSC_DETERMINE; 594f7320561SMichael Lange } 595c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cell faces %" PetscInt_FMT "\n", numCellFaces)); 596f7320561SMichael Lange } 59719d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 598f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 599930bae4bSMatthew G. Knepley numFaces = (PetscInt)(s.last - s.first + 1); 60035462f7fSMichael Lange if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 601895fcc3eSMichael Lange else numFaceVertices = s.nd; 602c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of faces %" PetscInt_FMT " face vertices: %" PetscInt_FMT "\n", numFaces, numFaceVertices)); 603f7320561SMichael Lange } else { /* Data section */ 604cf554e2cSMatthew G. Knepley unsigned int z; 605cf554e2cSMatthew G. Knepley 6061dca8a05SBarry Smith PetscCheck(numFaceVertices == PETSC_DETERMINE || s.nd == numFaceVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 60708401ef6SPierre Jolivet PetscCheck(numFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 608f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 609f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 6109566063dSJacob Faibussowitsch if (!faces) PetscCall(PetscMalloc1(numFaces * numFaceEntries, &faces)); 6119566063dSJacob Faibussowitsch if (!faceZoneIDs) PetscCall(PetscMalloc1(numFaces, &faceZoneIDs)); 6129566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&faces[(s.first - 1) * numFaceEntries], s.data, (s.last - s.first + 1) * numFaceEntries * sizeof(PetscInt))); 613ec78a56aSMichael Lange /* Record the zoneID for each face set */ 614cf554e2cSMatthew G. Knepley for (z = s.first - 1; z < s.last; z++) faceZoneIDs[z] = s.zoneID; 6159566063dSJacob Faibussowitsch PetscCall(PetscFree(s.data)); 616f7320561SMichael Lange } 617c7cc6ee6SMatthew G. Knepley } else if (s.index == 39) { /* Label information */ 618c7cc6ee6SMatthew G. Knepley if (s.zoneID >= maxZoneID) { 619c7cc6ee6SMatthew G. Knepley DMLabel *tmpL; 620c7cc6ee6SMatthew G. Knepley const char **tmp; 621c7cc6ee6SMatthew G. Knepley unsigned int newmax = maxZoneID + 1; 622c7cc6ee6SMatthew G. Knepley 623c7cc6ee6SMatthew G. Knepley while (newmax < s.zoneID + 1) newmax *= 2; 624c7cc6ee6SMatthew G. Knepley PetscCall(PetscCalloc2(newmax, &tmp, newmax, &tmpL)); 625c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) { 626c7cc6ee6SMatthew G. Knepley tmp[i] = zoneNames[i]; 627c7cc6ee6SMatthew G. Knepley tmpL[i] = zoneLabels[i]; 628c7cc6ee6SMatthew G. Knepley } 629c7cc6ee6SMatthew G. Knepley maxZoneID = newmax; 630c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 631c7cc6ee6SMatthew G. Knepley zoneNames = tmp; 632c7cc6ee6SMatthew G. Knepley zoneLabels = tmpL; 633c7cc6ee6SMatthew G. Knepley } 634c7cc6ee6SMatthew G. Knepley zoneNames[s.zoneID] = (const char *)s.data; 6352f0bd6dcSMichael Lange } 6362f0bd6dcSMichael Lange } while (s.index >= 0); 6372f0bd6dcSMichael Lange } 6389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm)); 63908401ef6SPierre Jolivet PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 640f7320561SMichael Lange 641f7320561SMichael Lange /* Allocate cell-vertex mesh */ 6429566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 6439566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 6449566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 645c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetChart(*dm, 0, numCells + numFaces + numVertices)); 646dd400576SPatrick Sanan if (rank == 0) { 64708401ef6SPierre Jolivet PetscCheck(numCells >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 648895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 649c7cc6ee6SMatthew G. Knepley if (numCellFaces == PETSC_DETERMINE) { 650c7cc6ee6SMatthew G. Knepley numCellFaces = numFaceVertices + 1; 651c7cc6ee6SMatthew G. Knepley ct = numCellFaces == 3 ? DM_POLYTOPE_TRIANGLE : (numCellFaces == 4 ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_UNKNOWN_CELL); 652c7cc6ee6SMatthew G. Knepley } 653c7cc6ee6SMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(*dm, c, numCellFaces)); 654c7cc6ee6SMatthew G. Knepley for (PetscInt f = 0; f < numFaces; ++f) PetscCall(DMPlexSetConeSize(*dm, f + numCells + numVertices, numFaceVertices)); 655f7320561SMichael Lange } 6569566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 657f7320561SMichael Lange 658dd400576SPatrick Sanan if (rank == 0 && faces) { 659c7cc6ee6SMatthew G. Knepley PetscInt *cones; 660f7320561SMichael Lange 661c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCones(*dm, &cones)); 662c7cc6ee6SMatthew G. Knepley PetscCheck(numCellFaces < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of cell faces %" PetscInt_FMT " exceeeds temporary storage", numCellFaces); 663c7cc6ee6SMatthew G. Knepley PetscCheck(numFaceVertices < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of face vertices %" PetscInt_FMT " exceeeds temporary storage", numFaceVertices); 664c7cc6ee6SMatthew G. Knepley for (PetscInt c = 0; c < numCells * numCellFaces; ++c) cones[c] = -1; 665c7cc6ee6SMatthew G. Knepley for (PetscInt f = 0; f < numFaces; f++) { 666c7cc6ee6SMatthew G. Knepley const PetscInt cl = faces[f * numFaceEntries + numFaceVertices] - 1; 667c7cc6ee6SMatthew G. Knepley const PetscInt cr = faces[f * numFaceEntries + numFaceVertices + 1] - 1; 668c7cc6ee6SMatthew G. Knepley const PetscInt face = f + numCells + numVertices; 669c7cc6ee6SMatthew G. Knepley const PetscInt *fc = &faces[f * numFaceEntries]; 670c7cc6ee6SMatthew G. Knepley PetscInt fcone[16]; 671c7cc6ee6SMatthew G. Knepley 672c7cc6ee6SMatthew G. Knepley // How could Fluent define the outward normal differently? Is there no end to the pain? 673c7cc6ee6SMatthew G. Knepley if (dim == 3) { 674c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, -1)); 675c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, 0)); 676c7cc6ee6SMatthew G. Knepley } else { 677c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, 0)); 678c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, -1)); 6799371c9d4SSatish Balay } 680c7cc6ee6SMatthew G. Knepley for (PetscInt v = 0; v < numFaceVertices; ++v) fcone[v] = fc[v] + numCells - 1; 681c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(*dm, face, fcone)); 682f7320561SMichael Lange } 683f7320561SMichael Lange } 6849566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 6859566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 686c7cc6ee6SMatthew G. Knepley if (dim == 3) { 6875fd9971aSMatthew G. Knepley DM idm; 688f7320561SMichael Lange 689c7cc6ee6SMatthew G. Knepley PetscCall(DMCreate(PetscObjectComm((PetscObject)*dm), &idm)); 690c7cc6ee6SMatthew G. Knepley PetscCall(DMSetType(idm, DMPLEX)); 691c7cc6ee6SMatthew G. Knepley PetscCall(DMSetDimension(idm, dim)); 692c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInterpolateFaces_Internal(*dm, 1, idm)); 6939566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 694f7320561SMichael Lange *dm = idm; 695f7320561SMichael Lange } 696c7cc6ee6SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-cas_dm_view")); 697c7cc6ee6SMatthew G. Knepley if (rank == 0 && faces) { 698*fc8f8d65SMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) PetscCall(ReorderCell(viewer, *dm, c, ct)); 699c7cc6ee6SMatthew G. Knepley } 700f7320561SMichael Lange 701dd400576SPatrick Sanan if (rank == 0 && faces) { 702631eb916SMichael Lange PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 703631eb916SMichael Lange const PetscInt *join, *meet; 7049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFaceVertices, &fverts)); 705ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 706ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 707631eb916SMichael Lange const PetscInt cl = faces[f * numFaceEntries + numFaceVertices] - 1; 708631eb916SMichael Lange const PetscInt cr = faces[f * numFaceEntries + numFaceVertices + 1] - 1; 709c7cc6ee6SMatthew G. Knepley const PetscInt id = (PetscInt)faceZoneIDs[f]; 710c7cc6ee6SMatthew G. Knepley 711631eb916SMichael Lange if (cl > 0 && cr > 0) { 712631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 7139371c9d4SSatish Balay cells[0] = cl; 7149371c9d4SSatish Balay cells[1] = cr; 7159566063dSJacob Faibussowitsch PetscCall(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet)); 716c7cc6ee6SMatthew 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); 717c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], id)); 718c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], meet[0], 1)); 7199566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet)); 720631eb916SMichael Lange } else { 721ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f * numFaceEntries + fi] + numCells - 1; 7229566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 72363a3b9bcSJacob Faibussowitsch PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f); 724c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], id)); 725c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], join[0], 1)); 7269566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 727ec78a56aSMichael Lange } 728631eb916SMichael Lange } 7299566063dSJacob Faibussowitsch PetscCall(PetscFree(fverts)); 730ec78a56aSMichael Lange } 731ec78a56aSMichael Lange 7327368db69SLisandro Dalcin { /* Create Face Sets label at all processes */ 7339371c9d4SSatish Balay enum { 7349371c9d4SSatish Balay n = 1 7359371c9d4SSatish Balay }; 7367368db69SLisandro Dalcin PetscBool flag[n]; 7377368db69SLisandro Dalcin 7387368db69SLisandro Dalcin flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE; 7399566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 7409566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 741c7cc6ee6SMatthew G. Knepley // TODO Code to create all the zone labels on each process 742c7cc6ee6SMatthew G. Knepley } 743c7cc6ee6SMatthew G. Knepley 744c7cc6ee6SMatthew G. Knepley if (!interpolate) { 745c7cc6ee6SMatthew G. Knepley DM udm; 746c7cc6ee6SMatthew G. Knepley 747c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexUninterpolate(*dm, &udm)); 748c7cc6ee6SMatthew G. Knepley PetscCall(DMDestroy(dm)); 749c7cc6ee6SMatthew G. Knepley *dm = udm; 7507368db69SLisandro Dalcin } 7517368db69SLisandro Dalcin 7523f6dc66eSMichael Lange /* Read coordinates */ 7539566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 7549566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 7559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 7569566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 757c7cc6ee6SMatthew G. Knepley for (PetscInt v = numCells; v < numCells + numVertices; ++v) { 7589566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dim)); 7599566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 7603f6dc66eSMichael Lange } 7619566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 7629566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 7639566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 7649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 7659566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 7669566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 7679566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 768dd400576SPatrick Sanan if (rank == 0 && coordsIn) { 769c7cc6ee6SMatthew G. Knepley for (PetscInt v = 0; v < numVertices; ++v) { 770c7cc6ee6SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) coords[v * dim + d] = coordsIn[v * dim + d]; 7713f6dc66eSMichael Lange } 7723f6dc66eSMichael Lange } 7739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 7749566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 7759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 7767368db69SLisandro Dalcin 777dd400576SPatrick Sanan if (rank == 0) { 7789566063dSJacob Faibussowitsch PetscCall(PetscFree(cellVertices)); 7799566063dSJacob Faibussowitsch PetscCall(PetscFree(faces)); 7809566063dSJacob Faibussowitsch PetscCall(PetscFree(faceZoneIDs)); 7819566063dSJacob Faibussowitsch PetscCall(PetscFree(coordsIn)); 782c7cc6ee6SMatthew G. Knepley if (zoneNames) 783c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) PetscCall(PetscFree(zoneNames[i])); 784c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 785f7320561SMichael Lange } 7863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7872f0bd6dcSMichael Lange } 788