12f0bd6dcSMichael Lange #define PETSCDM_DLL 22f0bd6dcSMichael Lange #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 32f0bd6dcSMichael Lange 42f0bd6dcSMichael Lange #undef __FUNCT__ 52f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluentFromFile" 62f0bd6dcSMichael Lange /*@C 72f0bd6dcSMichael Lange DMPlexCreateFluentFromFile - Create a DMPlex mesh from a Fluent mesh file 82f0bd6dcSMichael Lange 92f0bd6dcSMichael Lange + comm - The MPI communicator 102f0bd6dcSMichael Lange . filename - Name of the Fluent mesh file 112f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 122f0bd6dcSMichael Lange 132f0bd6dcSMichael Lange Output Parameter: 142f0bd6dcSMichael Lange . dm - The DM object representing the mesh 152f0bd6dcSMichael Lange 162f0bd6dcSMichael Lange Level: beginner 172f0bd6dcSMichael Lange 182f0bd6dcSMichael Lange .seealso: DMPlexCreateFromFile(), DMPlexCreateFluent(), DMPlexCreate() 192f0bd6dcSMichael Lange @*/ 202f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluentFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 212f0bd6dcSMichael Lange { 222f0bd6dcSMichael Lange PetscViewer viewer; 232f0bd6dcSMichael Lange PetscErrorCode ierr; 242f0bd6dcSMichael Lange 252f0bd6dcSMichael Lange PetscFunctionBegin; 262f0bd6dcSMichael Lange /* Create file viewer and build plex */ 272f0bd6dcSMichael Lange ierr = PetscViewerCreate(comm, &viewer);CHKERRQ(ierr); 282f0bd6dcSMichael Lange ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);CHKERRQ(ierr); 292f0bd6dcSMichael Lange ierr = PetscViewerFileSetMode(viewer, FILE_MODE_READ);CHKERRQ(ierr); 302f0bd6dcSMichael Lange ierr = PetscViewerFileSetName(viewer, filename);CHKERRQ(ierr); 312f0bd6dcSMichael Lange ierr = DMPlexCreateFluent(comm, viewer, interpolate, dm);CHKERRQ(ierr); 322f0bd6dcSMichael Lange ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 332f0bd6dcSMichael Lange PetscFunctionReturn(0); 342f0bd6dcSMichael Lange } 352f0bd6dcSMichael Lange 362f0bd6dcSMichael Lange #undef __FUNCT__ 372f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent_ReadString" 382f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent_ReadString(PetscViewer viewer, char *buffer, char delim) 392f0bd6dcSMichael Lange { 402f0bd6dcSMichael Lange PetscInt i = 0; 412f0bd6dcSMichael Lange PetscErrorCode ierr; 422f0bd6dcSMichael Lange 432f0bd6dcSMichael Lange PetscFunctionBegin; 442f0bd6dcSMichael Lange do {ierr = PetscViewerRead(viewer, &(buffer[i++]), 1, PETSC_CHAR);CHKERRQ(ierr);} 452f0bd6dcSMichael Lange while (buffer[i-1] != '\0' && buffer[i-1] != delim); 462f0bd6dcSMichael Lange buffer[i] = '\0'; 472f0bd6dcSMichael Lange PetscFunctionReturn(0); 482f0bd6dcSMichael Lange } 492f0bd6dcSMichael Lange 502f0bd6dcSMichael Lange #undef __FUNCT__ 512f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent_ReadSection" 522f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent_ReadSection(PetscViewer viewer, FluentSection *s) 532f0bd6dcSMichael Lange { 542f0bd6dcSMichael Lange char buffer[PETSC_MAX_PATH_LEN]; 552f0bd6dcSMichael Lange int snum; 562f0bd6dcSMichael Lange PetscErrorCode ierr; 572f0bd6dcSMichael Lange 582f0bd6dcSMichael Lange PetscFunctionBegin; 592f0bd6dcSMichael Lange /* Fast-forward to next section and derive its index */ 602f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 612f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ' ');CHKERRQ(ierr); 622f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->index)); 632f0bd6dcSMichael Lange /* If we can't match an index return -1 to signal end-of-file */ 642f0bd6dcSMichael Lange if (snum < 1) {s->index = -1; PetscFunctionReturn(0);} 652f0bd6dcSMichael Lange 662f0bd6dcSMichael Lange if (s->index == 0) { /* Comment */ 672f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 682f0bd6dcSMichael Lange 692f0bd6dcSMichael Lange } else if (s->index == 2) { /* Dimension */ 702f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 712f0bd6dcSMichael Lange snum = sscanf(buffer, "%d", &(s->nd)); 722f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 732f0bd6dcSMichael Lange 742f0bd6dcSMichael Lange } else if (s->index == 10) { /* Vertices */ 752f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 762f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 772f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 783f6dc66eSMichael Lange if (s->zoneID > 0) { 793f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 803f6dc66eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 813f6dc66eSMichael Lange ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); 823f6dc66eSMichael Lange ierr = PetscViewerRead(viewer, (PetscScalar*)s->data, numCoords*s->nd, PETSC_REAL);CHKERRQ(ierr); 833f6dc66eSMichael Lange } 842f0bd6dcSMichael Lange 852f0bd6dcSMichael Lange } else if (s->index == 12) { /* Cells */ 862f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 872f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 882f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 892f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 902f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 912f0bd6dcSMichael Lange if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 92f7320561SMichael Lange } else { /* Data section */ 932f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 942f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 95895fcc3eSMichael Lange if (s->nd == 0) { 96895fcc3eSMichael Lange /* Read cell type definitions for mixed cells */ 97895fcc3eSMichael Lange PetscInt i, numCells = s->last - s->first + 1; 98895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 99895fcc3eSMichael Lange for (i = 0; i < numCells; i++) { 100895fcc3eSMichael Lange ierr = PetscViewerRead(viewer, buffer, 1, PETSC_STRING);CHKERRQ(ierr); 101895fcc3eSMichael Lange } 102895fcc3eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 103895fcc3eSMichael Lange } 1042f0bd6dcSMichael Lange } 1052f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1062f0bd6dcSMichael Lange 1072f0bd6dcSMichael Lange } else if (s->index == 13) { /* Faces */ 1082f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1092f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1102f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1112f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1122f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1132f0bd6dcSMichael Lange if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 114f7320561SMichael Lange } else { /* Data section */ 115f7320561SMichael Lange PetscInt f, e, numEntries, numFaces; 116895fcc3eSMichael Lange int numFaceVert; 1172f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1182f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 119f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 120f7320561SMichael Lange switch (s->nd) { 121895fcc3eSMichael Lange case 0: numEntries = PETSC_DETERMINE; break; 122f7320561SMichael Lange case 2: numEntries = 2 + 2; break; /* linear */ 123f7320561SMichael Lange case 3: numEntries = 2 + 3; break; /* triangular */ 124f7320561SMichael Lange case 4: numEntries = 2 + 4; break; /* quadrilateral */ 125f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 126f7320561SMichael Lange } 127f7320561SMichael Lange numFaces = s->last-s->first + 1; 128895fcc3eSMichael Lange if (numEntries != PETSC_DETERMINE) { 129895fcc3eSMichael Lange /* Allocate space only if we already know the size of the block */ 130f7320561SMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 131895fcc3eSMichael Lange } 132f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 133895fcc3eSMichael Lange if (s->nd == 0) { 134895fcc3eSMichael Lange /* Determine the size of the block for "mixed" facets */ 135895fcc3eSMichael Lange ierr = PetscViewerRead(viewer, buffer, 1, PETSC_STRING);CHKERRQ(ierr); 136895fcc3eSMichael Lange snum = sscanf(buffer, "%x", &numFaceVert); 137895fcc3eSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 138895fcc3eSMichael Lange if (numEntries == PETSC_DETERMINE) { 139895fcc3eSMichael Lange numEntries = numFaceVert + 2; 140895fcc3eSMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 141895fcc3eSMichael Lange } else { 142895fcc3eSMichael Lange if (numEntries != numFaceVert + 2) { 143895fcc3eSMichael Lange SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for mixed faces in Fluent files"); 144895fcc3eSMichael Lange } 145895fcc3eSMichael Lange } 146895fcc3eSMichael Lange } 147f7320561SMichael Lange for (e = 0; e < numEntries; e++) { 148f7320561SMichael Lange ierr = PetscViewerRead(viewer, buffer, 1, PETSC_STRING);CHKERRQ(ierr); 149f7320561SMichael Lange snum = sscanf(buffer, "%x", &(((PetscInt*)s->data)[f*numEntries + e])); 150f7320561SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 151f7320561SMichael Lange } 152f7320561SMichael Lange } 153895fcc3eSMichael Lange s->nd = numEntries - 2; 154f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1552f0bd6dcSMichael Lange } 1562f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1572f0bd6dcSMichael Lange 1582f0bd6dcSMichael Lange } else { /* Unknown section type */ 1592f0bd6dcSMichael Lange PetscInt depth = 1; 1602f0bd6dcSMichael Lange do { 1612f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 1622f0bd6dcSMichael Lange do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, PETSC_CHAR);CHKERRQ(ierr);} 1632f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 1642f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 1652f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 1662f0bd6dcSMichael Lange } while (depth > 0); 1672f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); 1682f0bd6dcSMichael Lange } 1692f0bd6dcSMichael Lange PetscFunctionReturn(0); 1702f0bd6dcSMichael Lange } 1712f0bd6dcSMichael Lange 1722f0bd6dcSMichael Lange #undef __FUNCT__ 1732f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent" 1742f0bd6dcSMichael Lange /*@C 1752f0bd6dcSMichael Lange DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 1762f0bd6dcSMichael Lange 1772f0bd6dcSMichael Lange Collective on comm 1782f0bd6dcSMichael Lange 1792f0bd6dcSMichael Lange Input Parameters: 1802f0bd6dcSMichael Lange + comm - The MPI communicator 1812f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file 1822f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 1832f0bd6dcSMichael Lange 1842f0bd6dcSMichael Lange Output Parameter: 1852f0bd6dcSMichael Lange . dm - The DM object representing the mesh 1862f0bd6dcSMichael Lange 1872f0bd6dcSMichael Lange Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 1882f0bd6dcSMichael Lange 1892f0bd6dcSMichael Lange Level: beginner 1902f0bd6dcSMichael Lange 1912f0bd6dcSMichael Lange .keywords: mesh, fluent, case 1922f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate() 1932f0bd6dcSMichael Lange @*/ 1942f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1952f0bd6dcSMichael Lange { 1962f0bd6dcSMichael Lange PetscMPIInt rank; 197f7320561SMichael Lange PetscInt c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE; 198f7320561SMichael Lange PetscInt numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE; 199*ec78a56aSMichael Lange PetscInt *faces = NULL, *cellVertices, *faceZoneIDs = NULL; 2003f6dc66eSMichael Lange PetscInt d, coordSize; 2013f6dc66eSMichael Lange PetscScalar *coords, *coordsIn = NULL; 2023f6dc66eSMichael Lange PetscSection coordSection; 2033f6dc66eSMichael Lange Vec coordinates; 2042f0bd6dcSMichael Lange PetscErrorCode ierr; 2052f0bd6dcSMichael Lange 2062f0bd6dcSMichael Lange PetscFunctionBegin; 2072f0bd6dcSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 2082f0bd6dcSMichael Lange 2092f0bd6dcSMichael Lange if (!rank) { 2102f0bd6dcSMichael Lange FluentSection s; 2112f0bd6dcSMichael Lange do { 2122f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr); 2132f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 214f7320561SMichael Lange dim = s.nd; 2152f0bd6dcSMichael Lange 2162f0bd6dcSMichael Lange } else if (s.index == 10) { /* Vertices */ 217f7320561SMichael Lange if (s.zoneID == 0) numVertices = s.last; 2183f6dc66eSMichael Lange else { 2193f6dc66eSMichael Lange if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 2203f6dc66eSMichael Lange coordsIn = s.data; 2213f6dc66eSMichael Lange } 2222f0bd6dcSMichael Lange 2232f0bd6dcSMichael Lange } else if (s.index == 12) { /* Cells */ 224f7320561SMichael Lange if (s.zoneID == 0) numCells = s.last; 225f7320561SMichael Lange else { 226f7320561SMichael Lange switch (s.nd) { 227895fcc3eSMichael Lange case 0: numCellVertices = PETSC_DETERMINE; break; 228f7320561SMichael Lange case 1: numCellVertices = 3; break; /* triangular */ 229f7320561SMichael Lange case 2: numCellVertices = 4; break; /* tetrahedral */ 230f7320561SMichael Lange case 3: numCellVertices = 4; break; /* quadrilateral */ 231f7320561SMichael Lange case 4: numCellVertices = 8; break; /* hexahedral */ 232f7320561SMichael Lange case 5: numCellVertices = 5; break; /* pyramid */ 233f7320561SMichael Lange case 6: numCellVertices = 6; break; /* wedge */ 234895fcc3eSMichael Lange default: numCellVertices = PETSC_DETERMINE; 235f7320561SMichael Lange } 236f7320561SMichael Lange } 2372f0bd6dcSMichael Lange 2382f0bd6dcSMichael Lange } else if (s.index == 13) { /* Facets */ 239f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 240f7320561SMichael Lange numFaces = s.last - s.first + 1; 241895fcc3eSMichael Lange if (s.nd == 0) numFaceVertices = PETSC_DETERMINE; 242895fcc3eSMichael Lange else numFaceVertices = s.nd; 243f7320561SMichael Lange } else { /* Data section */ 244895fcc3eSMichael Lange if (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices) { 245895fcc3eSMichael Lange SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are not supported"); 246f7320561SMichael Lange } 247895fcc3eSMichael Lange if (numFaces < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 248f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 249f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 250f7320561SMichael Lange if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);} 251*ec78a56aSMichael Lange if (!faceZoneIDs) {ierr = PetscMalloc1(numFaces, &faceZoneIDs);CHKERRQ(ierr);} 252f7320561SMichael Lange ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr); 253*ec78a56aSMichael Lange /* Record the zoneID for each face set */ 254*ec78a56aSMichael Lange for (f = s.first -1; f < s.last; f++) faceZoneIDs[f] = s.zoneID; 255f7320561SMichael Lange ierr = PetscFree(s.data);CHKERRQ(ierr); 256f7320561SMichael Lange } 2572f0bd6dcSMichael Lange } 2582f0bd6dcSMichael Lange } while (s.index >= 0); 2592f0bd6dcSMichael Lange } 260f7320561SMichael Lange ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 261f7320561SMichael Lange if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 262f7320561SMichael Lange 263f7320561SMichael Lange /* Allocate cell-vertex mesh */ 264f7320561SMichael Lange ierr = DMCreate(comm, dm);CHKERRQ(ierr); 265f7320561SMichael Lange ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 266f7320561SMichael Lange ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 267f7320561SMichael Lange ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr); 268f7320561SMichael Lange if (!rank) { 269895fcc3eSMichael Lange if (numCells < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown number of cells in Fluent file"); 270895fcc3eSMichael Lange /* If no cell type was given we assume simplices */ 271895fcc3eSMichael Lange if (numCellVertices == PETSC_DETERMINE) numCellVertices = numFaceVertices + 1; 272f7320561SMichael Lange for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);} 273f7320561SMichael Lange } 274f7320561SMichael Lange ierr = DMSetUp(*dm);CHKERRQ(ierr); 275f7320561SMichael Lange 276f7320561SMichael Lange if (!rank) { 277f7320561SMichael Lange /* Derive cell-vertex list from face-vertex and face-cell maps */ 278f7320561SMichael Lange ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr); 279f7320561SMichael Lange for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 280f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 281f7320561SMichael Lange PetscInt *cell; 282f7320561SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 283f7320561SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 284f7320561SMichael Lange const PetscInt *face = &(faces[f*numFaceEntries]); 285f7320561SMichael Lange 286f7320561SMichael Lange if (cl > 0) { 287f7320561SMichael Lange cell = &(cellVertices[(cl-1) * numCellVertices]); 288f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 289f7320561SMichael Lange PetscBool found = PETSC_FALSE; 290f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 291f7320561SMichael Lange if (cell[c] < 0) break; 292f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 293f7320561SMichael Lange } 294f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 295f7320561SMichael Lange } 296f7320561SMichael Lange } 297f7320561SMichael Lange if (cr > 0) { 298f7320561SMichael Lange cell = &(cellVertices[(cr-1) * numCellVertices]); 299f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 300f7320561SMichael Lange PetscBool found = PETSC_FALSE; 301f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 302f7320561SMichael Lange if (cell[c] < 0) break; 303f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 304f7320561SMichael Lange } 305f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 306f7320561SMichael Lange } 307f7320561SMichael Lange } 308f7320561SMichael Lange } 309f7320561SMichael Lange for (c = 0; c < numCells; c++) { 310f7320561SMichael Lange ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr); 311f7320561SMichael Lange } 312f7320561SMichael Lange } 313f7320561SMichael Lange ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 314f7320561SMichael Lange ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 315f7320561SMichael Lange if (interpolate) { 316f7320561SMichael Lange DM idm = NULL; 317f7320561SMichael Lange 318f7320561SMichael Lange ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 319f7320561SMichael Lange ierr = DMDestroy(dm);CHKERRQ(ierr); 320f7320561SMichael Lange *dm = idm; 321f7320561SMichael Lange } 322f7320561SMichael Lange 323*ec78a56aSMichael Lange if (!rank) { 324*ec78a56aSMichael Lange PetscInt fi, *fverts; 325*ec78a56aSMichael Lange ierr = PetscMalloc1(numFaceVertices, &fverts);CHKERRQ(ierr); 326*ec78a56aSMichael Lange /* Mark facets by finding the full join of all adjacent vertices */ 327*ec78a56aSMichael Lange for (f = 0; f < numFaces; f++) { 328*ec78a56aSMichael Lange PetscInt joinSize; 329*ec78a56aSMichael Lange const PetscInt *join; 330*ec78a56aSMichael Lange for (fi = 0; fi < numFaceVertices; fi++) fverts[fi] = faces[f*numFaceEntries + fi] + numCells - 1; 331*ec78a56aSMichael Lange ierr = DMPlexGetFullJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 332*ec78a56aSMichael Lange if (joinSize != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not determine Plex facet for Fluent face %d", f); 333*ec78a56aSMichael Lange ierr = DMPlexSetLabelValue(*dm, "Face Sets", join[0], faceZoneIDs[f]);CHKERRQ(ierr); 334*ec78a56aSMichael Lange ierr = DMPlexRestoreJoin(*dm, numFaceVertices, fverts, &joinSize, &join);CHKERRQ(ierr); 335*ec78a56aSMichael Lange } 336*ec78a56aSMichael Lange ierr = PetscFree(fverts);CHKERRQ(ierr); 337*ec78a56aSMichael Lange } 338*ec78a56aSMichael Lange 3393f6dc66eSMichael Lange /* Read coordinates */ 3403f6dc66eSMichael Lange ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 3413f6dc66eSMichael Lange ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 3423f6dc66eSMichael Lange ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 3433f6dc66eSMichael Lange ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 3443f6dc66eSMichael Lange for (v = numCells; v < numCells+numVertices; ++v) { 3453f6dc66eSMichael Lange ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 3463f6dc66eSMichael Lange ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 3473f6dc66eSMichael Lange } 3483f6dc66eSMichael Lange ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 3493f6dc66eSMichael Lange ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 3503f6dc66eSMichael Lange ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 3513f6dc66eSMichael Lange ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 3523f6dc66eSMichael Lange ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 3533f6dc66eSMichael Lange ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 3543f6dc66eSMichael Lange ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 3553f6dc66eSMichael Lange if (!rank) { 3563f6dc66eSMichael Lange for (v = 0; v < numVertices; ++v) { 3573f6dc66eSMichael Lange for (d = 0; d < dim; ++d) { 3583f6dc66eSMichael Lange coords[v*dim+d] = coordsIn[v*dim+d]; 3593f6dc66eSMichael Lange } 3603f6dc66eSMichael Lange } 3613f6dc66eSMichael Lange } 3623f6dc66eSMichael Lange ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 3633f6dc66eSMichael Lange ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 3643f6dc66eSMichael Lange ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 365f7320561SMichael Lange if (!rank) { 366f7320561SMichael Lange ierr = PetscFree(cellVertices);CHKERRQ(ierr); 367f7320561SMichael Lange ierr = PetscFree(faces);CHKERRQ(ierr); 368*ec78a56aSMichael Lange ierr = PetscFree(faceZoneIDs);CHKERRQ(ierr); 3693f6dc66eSMichael Lange ierr = PetscFree(coordsIn);CHKERRQ(ierr); 370f7320561SMichael Lange } 3712f0bd6dcSMichael Lange PetscFunctionReturn(0); 3722f0bd6dcSMichael Lange } 373