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)); 43*c7cc6ee6SMatthew 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'; 46*c7cc6ee6SMatthew 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 50*c7cc6ee6SMatthew 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; 57*c7cc6ee6SMatthew 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++) { 71*c7cc6ee6SMatthew G. Knepley size_t len; 72*c7cc6ee6SMatthew 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; 77*c7cc6ee6SMatthew G. Knepley // Check for trailing parentheses 78*c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrlen(cbuf, &len)); 79*c7cc6ee6SMatthew G. Knepley while (cbuf[len - 1] == ')' && len > 0) { 80*c7cc6ee6SMatthew G. Knepley ++(*numClosingParens); 81*c7cc6ee6SMatthew G. Knepley --len; 82*c7cc6ee6SMatthew 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 */ 132*c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 133*c7cc6ee6SMatthew 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)); 141*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadValues(viewer, s->data, s->nd * numCoords, PETSC_SCALAR, s->index == 2010 ? PETSC_TRUE : PETSC_FALSE, &numClosingParens)); 142*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 143*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 1443f6dc66eSMichael Lange } 145*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 146*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 147*c7cc6ee6SMatthew 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 */ 149*c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 150*c7cc6ee6SMatthew 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)); 165*c7cc6ee6SMatthew 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)); 167*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 168*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 169895fcc3eSMichael Lange } 1702f0bd6dcSMichael Lange } 171*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 172*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 173*c7cc6ee6SMatthew 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 */ 175*c7cc6ee6SMatthew G. Knepley PetscInt numClosingParens = 0; 176*c7cc6ee6SMatthew 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; 213*c7cc6ee6SMatthew 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 } 221*c7cc6ee6SMatthew 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)); 224*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 225*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 2262f0bd6dcSMichael Lange } 227*c7cc6ee6SMatthew G. Knepley if (!numClosingParens) PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 228*c7cc6ee6SMatthew G. Knepley else --numClosingParens; 229*c7cc6ee6SMatthew G. Knepley PetscCheck(!numClosingParens, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 230*c7cc6ee6SMatthew G. Knepley } else if (s->index == 39) { /* Label information */ 231*c7cc6ee6SMatthew G. Knepley char labelName[PETSC_MAX_PATH_LEN]; 232*c7cc6ee6SMatthew G. Knepley char caseName[PETSC_MAX_PATH_LEN]; 2332f0bd6dcSMichael Lange 234*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, ')')); 235*c7cc6ee6SMatthew G. Knepley snum = sscanf(buffer, "(%d %s %s %d)", &s->zoneID, caseName, labelName, &s->nd); 236*c7cc6ee6SMatthew G. Knepley PetscCheck(snum == 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file: %d", snum); 237*c7cc6ee6SMatthew G. Knepley PetscInt depth = 1; 238*c7cc6ee6SMatthew G. Knepley do { 239*c7cc6ee6SMatthew G. Knepley /* Match parentheses when parsing unknown sections */ 240*c7cc6ee6SMatthew G. Knepley do PetscCall(PetscViewerRead(viewer, &buffer[0], 1, NULL, PETSC_CHAR)); 241*c7cc6ee6SMatthew G. Knepley while (buffer[0] != '(' && buffer[0] != ')'); 242*c7cc6ee6SMatthew G. Knepley if (buffer[0] == '(') depth++; 243*c7cc6ee6SMatthew G. Knepley if (buffer[0] == ')') depth--; 244*c7cc6ee6SMatthew G. Knepley } while (depth > 0); 245*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexCreateFluent_ReadString(viewer, buffer, '\n')); 246*c7cc6ee6SMatthew G. Knepley PetscCall(PetscStrallocpy(labelName, (char **)&s->data)); 247*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Zone ID %d 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 262*c7cc6ee6SMatthew G. Knepley // Inserts point `face` with orientation `ornt` into the cone of point `cell` at position `c`, which is the first empty slot 263*c7cc6ee6SMatthew G. Knepley static PetscErrorCode InsertFace(DM dm, PetscInt cell, PetscInt face, PetscInt ornt) 264*c7cc6ee6SMatthew G. Knepley { 265*c7cc6ee6SMatthew G. Knepley const PetscInt *cone; 266*c7cc6ee6SMatthew G. Knepley PetscInt coneSize, c; 267*c7cc6ee6SMatthew G. Knepley 268*c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 269*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cell, &cone)); 270*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 271*c7cc6ee6SMatthew G. Knepley for (c = 0; c < coneSize; ++c) 272*c7cc6ee6SMatthew G. Knepley if (cone[c] < 0) break; 273*c7cc6ee6SMatthew 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); 274*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertCone(dm, cell, c, face)); 275*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInsertConeOrientation(dm, cell, c, ornt)); 276*c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 277*c7cc6ee6SMatthew G. Knepley } 278*c7cc6ee6SMatthew G. Knepley 279*c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderPolygon(DM dm, PetscInt cell) 280*c7cc6ee6SMatthew G. Knepley { 281*c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt; 282*c7cc6ee6SMatthew G. Knepley PetscInt coneSize, newCone[16], newOrnt[16]; 283*c7cc6ee6SMatthew G. Knepley 284*c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 285*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 286*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dm, cell, &coneSize)); 287*c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 288*c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 289*c7cc6ee6SMatthew G. Knepley for (PetscInt c = 1; c < coneSize; ++c) { 290*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone; 291*c7cc6ee6SMatthew G. Knepley PetscInt firstVertex, lastVertex, c2; 292*c7cc6ee6SMatthew G. Knepley 293*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[c - 1], &fcone)); 294*c7cc6ee6SMatthew G. Knepley lastVertex = newOrnt[c - 1] ? fcone[0] : fcone[1]; 295*c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < coneSize; ++c2) { 296*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2; 297*c7cc6ee6SMatthew G. Knepley 298*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, cone[c2], &fcone2)); 299*c7cc6ee6SMatthew G. Knepley firstVertex = ornt[c2] ? fcone2[1] : fcone2[0]; 300*c7cc6ee6SMatthew G. Knepley if (lastVertex == firstVertex) { 301*c7cc6ee6SMatthew G. Knepley // Point `cell` matched point `lastVertex` on face `cone[c2]` with orientation `ornt[c2]` 302*c7cc6ee6SMatthew G. Knepley break; 303*c7cc6ee6SMatthew G. Knepley } 304*c7cc6ee6SMatthew G. Knepley } 305*c7cc6ee6SMatthew 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); 306*c7cc6ee6SMatthew G. Knepley newCone[c] = cone[c2]; 307*c7cc6ee6SMatthew G. Knepley newOrnt[c] = ornt[c2]; 308*c7cc6ee6SMatthew G. Knepley } 309*c7cc6ee6SMatthew G. Knepley { 310*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone, *fcone2; 311*c7cc6ee6SMatthew G. Knepley PetscInt vertex, vertex2; 312*c7cc6ee6SMatthew G. Knepley 313*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[coneSize - 1], &fcone)); 314*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCone(dm, newCone[0], &fcone2)); 315*c7cc6ee6SMatthew G. Knepley vertex = newOrnt[coneSize - 1] ? fcone[0] : fcone[1]; 316*c7cc6ee6SMatthew G. Knepley vertex2 = newOrnt[0] ? fcone2[1] : fcone2[0]; 317*c7cc6ee6SMatthew G. Knepley PetscCheck(vertex == vertex2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " did not match at the endpoint", cell); 318*c7cc6ee6SMatthew G. Knepley } 319*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 320*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 321*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 322*c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 323*c7cc6ee6SMatthew G. Knepley } 324*c7cc6ee6SMatthew G. Knepley 325*c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderTetrahedron(DM dm, PetscInt cell) 326*c7cc6ee6SMatthew G. Knepley { 327*c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr, faces[4] = {0, 1, 3, 2}; 328*c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 329*c7cc6ee6SMatthew G. Knepley 330*c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 331*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 332*c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 333*c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 334*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 335*c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, newOrnt[0]); 336*c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial triangle 337*c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 3; ++e) { 338*c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 339*c7cc6ee6SMatthew G. Knepley PetscInt c; 340*c7cc6ee6SMatthew G. Knepley 341*c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the tetrahedron 342*c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 343*c7cc6ee6SMatthew G. Knepley for (c = 1; c < 4; ++c) { 344*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 345*c7cc6ee6SMatthew G. Knepley PetscInt c2; 346*c7cc6ee6SMatthew G. Knepley 347*c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 348*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 349*c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_TRIANGLE, ornt[c]); 350*c7cc6ee6SMatthew G. Knepley // Check for edge 351*c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 3; ++c2) { 352*c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 353*c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 354*c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 355*c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 356*c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 357*c7cc6ee6SMatthew G. Knepley break; 358*c7cc6ee6SMatthew G. Knepley } 359*c7cc6ee6SMatthew G. Knepley } 360*c7cc6ee6SMatthew G. Knepley if (c2 < 3) { 361*c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 362*c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 363*c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_TRIANGLE, !e ? (c2 + 1) % 3 : c2, ornt[c]); 364*c7cc6ee6SMatthew G. Knepley break; 365*c7cc6ee6SMatthew G. Knepley } 366*c7cc6ee6SMatthew G. Knepley } 367*c7cc6ee6SMatthew 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); 368*c7cc6ee6SMatthew G. Knepley } 369*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 370*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 371*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 372*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 373*c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 374*c7cc6ee6SMatthew G. Knepley } 375*c7cc6ee6SMatthew G. Knepley 376*c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderHexahedron(DM dm, PetscInt cell) 377*c7cc6ee6SMatthew G. Knepley { 378*c7cc6ee6SMatthew G. Knepley const PetscInt *cone, *ornt, *fcone, *fornt, *farr; 379*c7cc6ee6SMatthew G. Knepley const PetscInt faces[6] = {0, 5, 3, 4, 2, 1}; 380*c7cc6ee6SMatthew G. Knepley PetscInt used[6] = {1, 0, 0, 0, 0, 0}; 381*c7cc6ee6SMatthew G. Knepley PetscInt newCone[16], newOrnt[16]; 382*c7cc6ee6SMatthew G. Knepley 383*c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 384*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cell, &cone, &ornt)); 385*c7cc6ee6SMatthew G. Knepley newCone[0] = cone[0]; 386*c7cc6ee6SMatthew G. Knepley newOrnt[0] = ornt[0]; 387*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[0], &fcone, &fornt)); 388*c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[0]); 389*c7cc6ee6SMatthew G. Knepley // Loop over each edge in the initial quadrilateral 390*c7cc6ee6SMatthew G. Knepley for (PetscInt e = 0; e < 4; ++e) { 391*c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 392*c7cc6ee6SMatthew G. Knepley PetscInt c; 393*c7cc6ee6SMatthew G. Knepley 394*c7cc6ee6SMatthew G. Knepley // Loop over each remaining face in the hexahedron 395*c7cc6ee6SMatthew G. Knepley // On face `newCone[0]`, trying to match edge `edge` with final orientation `eornt` to an edge on another face 396*c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) { 397*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 398*c7cc6ee6SMatthew G. Knepley PetscInt c2; 399*c7cc6ee6SMatthew G. Knepley 400*c7cc6ee6SMatthew G. Knepley // Checking face `cone[c]` with orientation `ornt[c]` 401*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 402*c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 403*c7cc6ee6SMatthew G. Knepley // Check for edge 404*c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 405*c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 406*c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 407*c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 408*c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 409*c7cc6ee6SMatthew G. Knepley // Matched face `newCone[0]` with orientation `newOrnt[0]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 410*c7cc6ee6SMatthew G. Knepley break; 411*c7cc6ee6SMatthew G. Knepley } 412*c7cc6ee6SMatthew G. Knepley } 413*c7cc6ee6SMatthew G. Knepley if (c2 < 4) { 414*c7cc6ee6SMatthew G. Knepley used[c] = 1; 415*c7cc6ee6SMatthew G. Knepley newCone[faces[e + 1]] = cone[c]; 416*c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched (only the first edge matches a side different from 0) 417*c7cc6ee6SMatthew G. Knepley newOrnt[faces[e + 1]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, !e ? (c2 + 1) % 4 : c2, ornt[c]); 418*c7cc6ee6SMatthew G. Knepley break; 419*c7cc6ee6SMatthew G. Knepley } 420*c7cc6ee6SMatthew G. Knepley } 421*c7cc6ee6SMatthew 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); 422*c7cc6ee6SMatthew G. Knepley } 423*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 424*c7cc6ee6SMatthew G. Knepley // Add last face 425*c7cc6ee6SMatthew G. Knepley { 426*c7cc6ee6SMatthew G. Knepley PetscInt c, c2; 427*c7cc6ee6SMatthew G. Knepley 428*c7cc6ee6SMatthew G. Knepley for (c = 1; c < 6; ++c) 429*c7cc6ee6SMatthew G. Knepley if (!used[c]) break; 430*c7cc6ee6SMatthew G. Knepley PetscCheck(c < 6, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell %" PetscInt_FMT " could not find an available face", cell); 431*c7cc6ee6SMatthew G. Knepley // Match first edge to 3rd edge in newCone[2] 432*c7cc6ee6SMatthew G. Knepley { 433*c7cc6ee6SMatthew G. Knepley const PetscInt *fcone2, *fornt2, *farr2; 434*c7cc6ee6SMatthew G. Knepley 435*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, newCone[2], &fcone, &fornt)); 436*c7cc6ee6SMatthew G. Knepley farr = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, newOrnt[2]); 437*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dm, cone[c], &fcone2, &fornt2)); 438*c7cc6ee6SMatthew G. Knepley farr2 = DMPolytopeTypeGetArrangement(DM_POLYTOPE_QUADRILATERAL, ornt[c]); 439*c7cc6ee6SMatthew G. Knepley 440*c7cc6ee6SMatthew G. Knepley const PetscInt e = 2; 441*c7cc6ee6SMatthew G. Knepley const PetscInt edge = fcone[farr[e * 2 + 0]], eornt = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr[e * 2 + 1], fornt[farr[e * 2 + 0]]); 442*c7cc6ee6SMatthew 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]` 443*c7cc6ee6SMatthew G. Knepley for (c2 = 0; c2 < 4; ++c2) { 444*c7cc6ee6SMatthew G. Knepley const PetscInt edge2 = fcone2[farr2[c2 * 2 + 0]], eornt2 = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_SEGMENT, farr2[c2 * 2 + 1], fornt2[farr2[c2 * 2 + 0]]); 445*c7cc6ee6SMatthew G. Knepley // Trying to match edge `edge2` with final orientation `eornt2` 446*c7cc6ee6SMatthew G. Knepley if (edge == edge2) { 447*c7cc6ee6SMatthew G. Knepley PetscCheck(eornt == -(eornt2 + 1), PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Edge %" PetscInt_FMT " found twice with the same orientation", edge); 448*c7cc6ee6SMatthew G. Knepley // Matched face `newCone[2]` with orientation `newOrnt[2]` to face `cone[c]` with orientation `ornt[c]` along edge `edge` 449*c7cc6ee6SMatthew G. Knepley break; 450*c7cc6ee6SMatthew G. Knepley } 451*c7cc6ee6SMatthew G. Knepley } 452*c7cc6ee6SMatthew G. Knepley PetscCheck(c2 < 4, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not fit last face in"); 453*c7cc6ee6SMatthew G. Knepley } 454*c7cc6ee6SMatthew G. Knepley newCone[faces[5]] = cone[c]; 455*c7cc6ee6SMatthew G. Knepley // Compute new orientation of face based on which edge was matched 456*c7cc6ee6SMatthew G. Knepley newOrnt[faces[5]] = DMPolytopeTypeComposeOrientation(DM_POLYTOPE_QUADRILATERAL, c2, ornt[c]); 457*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, newCone[0], &fcone, &fornt)); 458*c7cc6ee6SMatthew G. Knepley } 459*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(dm, cell, newCone)); 460*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetConeOrientation(dm, cell, newOrnt)); 461*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dm, cell, &cone, &ornt)); 462*c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 463*c7cc6ee6SMatthew G. Knepley } 464*c7cc6ee6SMatthew G. Knepley 465*c7cc6ee6SMatthew G. Knepley static PetscErrorCode ReorderCell(DM dm, PetscInt cell, DMPolytopeType ct) 466*c7cc6ee6SMatthew G. Knepley { 467*c7cc6ee6SMatthew G. Knepley PetscFunctionBegin; 468*c7cc6ee6SMatthew G. Knepley switch (ct) { 469*c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TRIANGLE: 470*c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_QUADRILATERAL: 471*c7cc6ee6SMatthew G. Knepley PetscCall(ReorderPolygon(dm, cell)); 472*c7cc6ee6SMatthew G. Knepley break; 473*c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_TETRAHEDRON: 474*c7cc6ee6SMatthew G. Knepley PetscCall(ReorderTetrahedron(dm, cell)); 475*c7cc6ee6SMatthew G. Knepley break; 476*c7cc6ee6SMatthew G. Knepley case DM_POLYTOPE_HEXAHEDRON: 477*c7cc6ee6SMatthew G. Knepley PetscCall(ReorderHexahedron(dm, cell)); 478*c7cc6ee6SMatthew G. Knepley break; 479*c7cc6ee6SMatthew G. Knepley default: 480*c7cc6ee6SMatthew G. Knepley PetscCheck(0, PETSC_COMM_SELF, PETSC_ERR_SUP, "Celltype %s is unsupported", DMPolytopeTypes[ct]); 481*c7cc6ee6SMatthew G. Knepley break; 482*c7cc6ee6SMatthew G. Knepley } 483*c7cc6ee6SMatthew G. Knepley PetscFunctionReturn(PETSC_SUCCESS); 484*c7cc6ee6SMatthew G. Knepley } 485*c7cc6ee6SMatthew G. Knepley 4862f0bd6dcSMichael Lange /*@C 4871d27aa22SBarry Smith DMPlexCreateFluent - Create a `DMPLEX` mesh from a Fluent mesh file <http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm>. 4882f0bd6dcSMichael Lange 489d083f849SBarry Smith Collective 4902f0bd6dcSMichael Lange 4912f0bd6dcSMichael Lange Input Parameters: 4922f0bd6dcSMichael Lange + comm - The MPI communicator 493a1cb98faSBarry Smith . viewer - The `PetscViewer` associated with a Fluent mesh file 4942f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 4952f0bd6dcSMichael Lange 4962f0bd6dcSMichael Lange Output Parameter: 497a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 4982f0bd6dcSMichael Lange 4992f0bd6dcSMichael Lange Level: beginner 5002f0bd6dcSMichael Lange 5011cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreate()` 5022f0bd6dcSMichael Lange @*/ 503d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 504d71ae5a4SJacob Faibussowitsch { 505*c7cc6ee6SMatthew G. Knepley PetscInt dim = PETSC_DETERMINE; 506*c7cc6ee6SMatthew G. Knepley PetscInt numCells = 0; 507*c7cc6ee6SMatthew G. Knepley PetscInt numVertices = 0; 508*c7cc6ee6SMatthew G. Knepley PetscInt numCellFaces = PETSC_DETERMINE; 509*c7cc6ee6SMatthew G. Knepley DMPolytopeType ct = DM_POLYTOPE_UNKNOWN_CELL; 510*c7cc6ee6SMatthew G. Knepley PetscInt numFaces = 0; 511*c7cc6ee6SMatthew G. Knepley PetscInt numFaceEntries = PETSC_DETERMINE; 512*c7cc6ee6SMatthew G. Knepley PetscInt numFaceVertices = PETSC_DETERMINE; 513*c7cc6ee6SMatthew G. Knepley PetscInt *faces = NULL; 514*c7cc6ee6SMatthew G. Knepley PetscInt *cellVertices = NULL; 515*c7cc6ee6SMatthew G. Knepley unsigned int *faceZoneIDs = NULL; 5167368db69SLisandro Dalcin DMLabel faceSets = NULL; 517*c7cc6ee6SMatthew G. Knepley DMLabel *zoneLabels = NULL; 518*c7cc6ee6SMatthew G. Knepley const char **zoneNames = NULL; 519*c7cc6ee6SMatthew G. Knepley unsigned int maxZoneID = 0; 520*c7cc6ee6SMatthew G. Knepley PetscScalar *coordsIn = NULL; 521*c7cc6ee6SMatthew G. Knepley PetscScalar *coords; 5223f6dc66eSMichael Lange PetscSection coordSection; 5233f6dc66eSMichael Lange Vec coordinates; 524*c7cc6ee6SMatthew G. Knepley PetscInt coordSize, f; 525*c7cc6ee6SMatthew G. Knepley PetscMPIInt rank; 5262f0bd6dcSMichael Lange 5272f0bd6dcSMichael Lange PetscFunctionBegin; 5289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5292f0bd6dcSMichael Lange 530dd400576SPatrick Sanan if (rank == 0) { 5312f0bd6dcSMichael Lange FluentSection s; 532*c7cc6ee6SMatthew G. Knepley numFaces = PETSC_DETERMINE; 5332f0bd6dcSMichael Lange do { 5349566063dSJacob Faibussowitsch PetscCall(DMPlexCreateFluent_ReadSection(viewer, &s)); 5352f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 536f7320561SMichael Lange dim = s.nd; 537*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found dimension: %" PetscInt_FMT "\n", dim)); 53819d58f9dSMichael Lange } else if (s.index == 10 || s.index == 2010) { /* Vertices */ 539*c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 540*c7cc6ee6SMatthew G. Knepley numVertices = s.last; 541*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of vertices: %" PetscInt_FMT "\n", numVertices)); 542*c7cc6ee6SMatthew G. Knepley } else { 543*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found vertex coordinates\n")); 54428b400f6SJacob Faibussowitsch PetscCheck(!coordsIn, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 5451a9c30ecSMatthew G. Knepley coordsIn = (PetscScalar *)s.data; 5463f6dc66eSMichael Lange } 5472f0bd6dcSMichael Lange 54819d58f9dSMichael Lange } else if (s.index == 12 || s.index == 2012) { /* Cells */ 549*c7cc6ee6SMatthew G. Knepley if (s.zoneID == 0) { 550*c7cc6ee6SMatthew G. Knepley numCells = s.last; 551*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cells %" PetscInt_FMT "\n", numCells)); 552*c7cc6ee6SMatthew G. Knepley } else { 553f7320561SMichael Lange switch (s.nd) { 554d71ae5a4SJacob Faibussowitsch case 0: 555*c7cc6ee6SMatthew G. Knepley numCellFaces = PETSC_DETERMINE; 556*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_POINT; 557d71ae5a4SJacob Faibussowitsch break; 558d71ae5a4SJacob Faibussowitsch case 1: 559*c7cc6ee6SMatthew G. Knepley numCellFaces = 3; 560*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TRIANGLE; 561*c7cc6ee6SMatthew G. Knepley break; 562d71ae5a4SJacob Faibussowitsch case 2: 563*c7cc6ee6SMatthew G. Knepley numCellFaces = 4; 564*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TETRAHEDRON; 565*c7cc6ee6SMatthew G. Knepley break; 566d71ae5a4SJacob Faibussowitsch case 3: 567*c7cc6ee6SMatthew G. Knepley numCellFaces = 4; 568*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_QUADRILATERAL; 569*c7cc6ee6SMatthew G. Knepley break; 570d71ae5a4SJacob Faibussowitsch case 4: 571*c7cc6ee6SMatthew G. Knepley numCellFaces = 6; 572*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_HEXAHEDRON; 573*c7cc6ee6SMatthew G. Knepley break; 574d71ae5a4SJacob Faibussowitsch case 5: 575*c7cc6ee6SMatthew G. Knepley numCellFaces = 5; 576*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_PYRAMID; 577*c7cc6ee6SMatthew G. Knepley break; 578d71ae5a4SJacob Faibussowitsch case 6: 579*c7cc6ee6SMatthew G. Knepley numCellFaces = 5; 580*c7cc6ee6SMatthew G. Knepley ct = DM_POLYTOPE_TRI_PRISM; 581*c7cc6ee6SMatthew G. Knepley break; 582d71ae5a4SJacob Faibussowitsch default: 583*c7cc6ee6SMatthew G. Knepley numCellFaces = PETSC_DETERMINE; 584f7320561SMichael Lange } 585*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of cell faces %" PetscInt_FMT "\n", numCellFaces)); 586f7320561SMichael Lange } 58719d58f9dSMichael Lange } else if (s.index == 13 || s.index == 2013) { /* Facets */ 588f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 589930bae4bSMatthew G. Knepley numFaces = (PetscInt)(s.last - s.first + 1); 59035462f7fSMichael Lange if (s.nd == 0 || s.nd == 5) numFaceVertices = PETSC_DETERMINE; 591895fcc3eSMichael Lange else numFaceVertices = s.nd; 592*c7cc6ee6SMatthew G. Knepley PetscCall(PetscInfo((PetscObject)viewer, "CASE: Found number of faces %" PetscInt_FMT " face vertices: %" PetscInt_FMT "\n", numFaces, numFaceVertices)); 593f7320561SMichael Lange } else { /* Data section */ 594cf554e2cSMatthew G. Knepley unsigned int z; 595cf554e2cSMatthew G. Knepley 5961dca8a05SBarry Smith PetscCheck(numFaceVertices == PETSC_DETERMINE || s.nd == numFaceVertices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 59708401ef6SPierre Jolivet PetscCheck(numFaces >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 598f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 599f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 6009566063dSJacob Faibussowitsch if (!faces) PetscCall(PetscMalloc1(numFaces * numFaceEntries, &faces)); 6019566063dSJacob Faibussowitsch if (!faceZoneIDs) PetscCall(PetscMalloc1(numFaces, &faceZoneIDs)); 6029566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(&faces[(s.first - 1) * numFaceEntries], s.data, (s.last - s.first + 1) * numFaceEntries * sizeof(PetscInt))); 603ec78a56aSMichael Lange /* Record the zoneID for each face set */ 604cf554e2cSMatthew G. Knepley for (z = s.first - 1; z < s.last; z++) faceZoneIDs[z] = s.zoneID; 6059566063dSJacob Faibussowitsch PetscCall(PetscFree(s.data)); 606f7320561SMichael Lange } 607*c7cc6ee6SMatthew G. Knepley } else if (s.index == 39) { /* Label information */ 608*c7cc6ee6SMatthew G. Knepley if (s.zoneID >= maxZoneID) { 609*c7cc6ee6SMatthew G. Knepley DMLabel *tmpL; 610*c7cc6ee6SMatthew G. Knepley const char **tmp; 611*c7cc6ee6SMatthew G. Knepley unsigned int newmax = maxZoneID + 1; 612*c7cc6ee6SMatthew G. Knepley 613*c7cc6ee6SMatthew G. Knepley while (newmax < s.zoneID + 1) newmax *= 2; 614*c7cc6ee6SMatthew G. Knepley PetscCall(PetscCalloc2(newmax, &tmp, newmax, &tmpL)); 615*c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) { 616*c7cc6ee6SMatthew G. Knepley tmp[i] = zoneNames[i]; 617*c7cc6ee6SMatthew G. Knepley tmpL[i] = zoneLabels[i]; 618*c7cc6ee6SMatthew G. Knepley } 619*c7cc6ee6SMatthew G. Knepley maxZoneID = newmax; 620*c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 621*c7cc6ee6SMatthew G. Knepley zoneNames = tmp; 622*c7cc6ee6SMatthew G. Knepley zoneLabels = tmpL; 623*c7cc6ee6SMatthew G. Knepley } 624*c7cc6ee6SMatthew G. Knepley zoneNames[s.zoneID] = (const char *)s.data; 6252f0bd6dcSMichael Lange } 6262f0bd6dcSMichael Lange } while (s.index >= 0); 6272f0bd6dcSMichael Lange } 6289566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPIU_INT, 0, comm)); 62908401ef6SPierre Jolivet PetscCheck(dim >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 630f7320561SMichael Lange 631f7320561SMichael Lange /* Allocate cell-vertex mesh */ 6329566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 6339566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 6349566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, dim)); 635*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetChart(*dm, 0, numCells + numFaces + numVertices)); 636dd400576SPatrick Sanan if (rank == 0) { 63708401ef6SPierre Jolivet PetscCheck(numCells >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 638895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 639*c7cc6ee6SMatthew G. Knepley if (numCellFaces == PETSC_DETERMINE) { 640*c7cc6ee6SMatthew G. Knepley numCellFaces = numFaceVertices + 1; 641*c7cc6ee6SMatthew G. Knepley ct = numCellFaces == 3 ? DM_POLYTOPE_TRIANGLE : (numCellFaces == 4 ? DM_POLYTOPE_TETRAHEDRON : DM_POLYTOPE_UNKNOWN_CELL); 642*c7cc6ee6SMatthew G. Knepley } 643*c7cc6ee6SMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) PetscCall(DMPlexSetConeSize(*dm, c, numCellFaces)); 644*c7cc6ee6SMatthew G. Knepley for (PetscInt f = 0; f < numFaces; ++f) PetscCall(DMPlexSetConeSize(*dm, f + numCells + numVertices, numFaceVertices)); 645f7320561SMichael Lange } 6469566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 647f7320561SMichael Lange 648dd400576SPatrick Sanan if (rank == 0 && faces) { 649*c7cc6ee6SMatthew G. Knepley PetscInt *cones; 650f7320561SMichael Lange 651*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexGetCones(*dm, &cones)); 652*c7cc6ee6SMatthew G. Knepley PetscCheck(numCellFaces < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of cell faces %" PetscInt_FMT " exceeeds temporary storage", numCellFaces); 653*c7cc6ee6SMatthew G. Knepley PetscCheck(numFaceVertices < 16, PETSC_COMM_SELF, PETSC_ERR_SUP, "Number of face vertices %" PetscInt_FMT " exceeeds temporary storage", numFaceVertices); 654*c7cc6ee6SMatthew G. Knepley for (PetscInt c = 0; c < numCells * numCellFaces; ++c) cones[c] = -1; 655*c7cc6ee6SMatthew G. Knepley for (PetscInt f = 0; f < numFaces; f++) { 656*c7cc6ee6SMatthew G. Knepley const PetscInt cl = faces[f * numFaceEntries + numFaceVertices] - 1; 657*c7cc6ee6SMatthew G. Knepley const PetscInt cr = faces[f * numFaceEntries + numFaceVertices + 1] - 1; 658*c7cc6ee6SMatthew G. Knepley const PetscInt face = f + numCells + numVertices; 659*c7cc6ee6SMatthew G. Knepley const PetscInt *fc = &faces[f * numFaceEntries]; 660*c7cc6ee6SMatthew G. Knepley PetscInt fcone[16]; 661*c7cc6ee6SMatthew G. Knepley 662*c7cc6ee6SMatthew G. Knepley // How could Fluent define the outward normal differently? Is there no end to the pain? 663*c7cc6ee6SMatthew G. Knepley if (dim == 3) { 664*c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, -1)); 665*c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, 0)); 666*c7cc6ee6SMatthew G. Knepley } else { 667*c7cc6ee6SMatthew G. Knepley if (cl >= 0) PetscCall(InsertFace(*dm, cl, face, 0)); 668*c7cc6ee6SMatthew G. Knepley if (cr >= 0) PetscCall(InsertFace(*dm, cr, face, -1)); 6699371c9d4SSatish Balay } 670*c7cc6ee6SMatthew G. Knepley for (PetscInt v = 0; v < numFaceVertices; ++v) fcone[v] = fc[v] + numCells - 1; 671*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexSetCone(*dm, face, fcone)); 672f7320561SMichael Lange } 673f7320561SMichael Lange } 6749566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 6759566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 676*c7cc6ee6SMatthew G. Knepley if (dim == 3) { 6775fd9971aSMatthew G. Knepley DM idm; 678f7320561SMichael Lange 679*c7cc6ee6SMatthew G. Knepley PetscCall(DMCreate(PetscObjectComm((PetscObject)*dm), &idm)); 680*c7cc6ee6SMatthew G. Knepley PetscCall(DMSetType(idm, DMPLEX)); 681*c7cc6ee6SMatthew G. Knepley PetscCall(DMSetDimension(idm, dim)); 682*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexInterpolateFaces_Internal(*dm, 1, idm)); 6839566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 684f7320561SMichael Lange *dm = idm; 685f7320561SMichael Lange } 686*c7cc6ee6SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-cas_dm_view")); 687*c7cc6ee6SMatthew G. Knepley if (rank == 0 && faces) { 688*c7cc6ee6SMatthew G. Knepley for (PetscInt c = 0; c < numCells; ++c) PetscCall(ReorderCell(*dm, c, ct)); 689*c7cc6ee6SMatthew G. Knepley } 690f7320561SMichael Lange 691dd400576SPatrick Sanan if (rank == 0 && faces) { 692631eb916SMichael Lange PetscInt fi, joinSize, meetSize, *fverts, cells[2]; 693631eb916SMichael Lange const PetscInt *join, *meet; 6949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numFaceVertices, &fverts)); 695ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 696ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 697631eb916SMichael Lange const PetscInt cl = faces[f * numFaceEntries + numFaceVertices] - 1; 698631eb916SMichael Lange const PetscInt cr = faces[f * numFaceEntries + numFaceVertices + 1] - 1; 699*c7cc6ee6SMatthew G. Knepley const PetscInt id = (PetscInt)faceZoneIDs[f]; 700*c7cc6ee6SMatthew G. Knepley 701631eb916SMichael Lange if (cl > 0 && cr > 0) { 702631eb916SMichael Lange /* If we know both adjoining cells we can use a single-level meet */ 7039371c9d4SSatish Balay cells[0] = cl; 7049371c9d4SSatish Balay cells[1] = cr; 7059566063dSJacob Faibussowitsch PetscCall(DMPlexGetMeet(*dm, 2, cells, &meetSize, &meet)); 706*c7cc6ee6SMatthew 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); 707*c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", meet[0], id)); 708*c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], meet[0], 1)); 7099566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreMeet(*dm, numFaceVertices, fverts, &meetSize, &meet)); 710631eb916SMichael Lange } else { 711ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f * numFaceEntries + fi] + numCells - 1; 7129566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 71363a3b9bcSJacob Faibussowitsch PetscCheck(joinSize == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %" PetscInt_FMT, f); 714*c7cc6ee6SMatthew G. Knepley PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", join[0], id)); 715*c7cc6ee6SMatthew G. Knepley if (zoneNames && zoneNames[id]) PetscCall(DMSetLabelValue_Fast(*dm, &zoneLabels[id], zoneNames[id], join[0], 1)); 7169566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join)); 717ec78a56aSMichael Lange } 718631eb916SMichael Lange } 7199566063dSJacob Faibussowitsch PetscCall(PetscFree(fverts)); 720ec78a56aSMichael Lange } 721ec78a56aSMichael Lange 7227368db69SLisandro Dalcin { /* Create Face Sets label at all processes */ 7239371c9d4SSatish Balay enum { 7249371c9d4SSatish Balay n = 1 7259371c9d4SSatish Balay }; 7267368db69SLisandro Dalcin PetscBool flag[n]; 7277368db69SLisandro Dalcin 7287368db69SLisandro Dalcin flag[0] = faceSets ? PETSC_TRUE : PETSC_FALSE; 7299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 7309566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 731*c7cc6ee6SMatthew G. Knepley // TODO Code to create all the zone labels on each process 732*c7cc6ee6SMatthew G. Knepley } 733*c7cc6ee6SMatthew G. Knepley 734*c7cc6ee6SMatthew G. Knepley if (!interpolate) { 735*c7cc6ee6SMatthew G. Knepley DM udm; 736*c7cc6ee6SMatthew G. Knepley 737*c7cc6ee6SMatthew G. Knepley PetscCall(DMPlexUninterpolate(*dm, &udm)); 738*c7cc6ee6SMatthew G. Knepley PetscCall(DMDestroy(dm)); 739*c7cc6ee6SMatthew G. Knepley *dm = udm; 7407368db69SLisandro Dalcin } 7417368db69SLisandro Dalcin 7423f6dc66eSMichael Lange /* Read coordinates */ 7439566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 7449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 7459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dim)); 7469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 747*c7cc6ee6SMatthew G. Knepley for (PetscInt v = numCells; v < numCells + numVertices; ++v) { 7489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dim)); 7499566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dim)); 7503f6dc66eSMichael Lange } 7519566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 7529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 7539566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 7549566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 7559566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 7569566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 7579566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 758dd400576SPatrick Sanan if (rank == 0 && coordsIn) { 759*c7cc6ee6SMatthew G. Knepley for (PetscInt v = 0; v < numVertices; ++v) { 760*c7cc6ee6SMatthew G. Knepley for (PetscInt d = 0; d < dim; ++d) coords[v * dim + d] = coordsIn[v * dim + d]; 7613f6dc66eSMichael Lange } 7623f6dc66eSMichael Lange } 7639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 7649566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 7659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 7667368db69SLisandro Dalcin 767dd400576SPatrick Sanan if (rank == 0) { 7689566063dSJacob Faibussowitsch PetscCall(PetscFree(cellVertices)); 7699566063dSJacob Faibussowitsch PetscCall(PetscFree(faces)); 7709566063dSJacob Faibussowitsch PetscCall(PetscFree(faceZoneIDs)); 7719566063dSJacob Faibussowitsch PetscCall(PetscFree(coordsIn)); 772*c7cc6ee6SMatthew G. Knepley if (zoneNames) 773*c7cc6ee6SMatthew G. Knepley for (PetscInt i = 0; i < (PetscInt)maxZoneID; ++i) PetscCall(PetscFree(zoneNames[i])); 774*c7cc6ee6SMatthew G. Knepley PetscCall(PetscFree2(zoneNames, zoneLabels)); 775f7320561SMichael Lange } 7763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7772f0bd6dcSMichael Lange } 778