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"); 78*3f6dc66eSMichael Lange if (s->zoneID > 0) { 79*3f6dc66eSMichael Lange PetscInt numCoords = s->last - s->first + 1; 80*3f6dc66eSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 81*3f6dc66eSMichael Lange ierr = PetscMalloc1(s->nd*numCoords, (PetscScalar**)&s->data);CHKERRQ(ierr); 82*3f6dc66eSMichael Lange ierr = PetscViewerRead(viewer, (PetscScalar*)s->data, numCoords*s->nd, PETSC_REAL);CHKERRQ(ierr); 83*3f6dc66eSMichael 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"); 952f0bd6dcSMichael Lange } 962f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 972f0bd6dcSMichael Lange 982f0bd6dcSMichael Lange } else if (s->index == 13) { /* Faces */ 992f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1002f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x", &(s->zoneID)); 1012f0bd6dcSMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 1022f0bd6dcSMichael Lange if (s->zoneID == 0) { /* Header section */ 1032f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd)); 1042f0bd6dcSMichael Lange if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 105f7320561SMichael Lange } else { /* Data section */ 106f7320561SMichael Lange PetscInt f, e, numEntries, numFaces; 1072f0bd6dcSMichael Lange snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd)); 1082f0bd6dcSMichael Lange if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 109f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr); 110f7320561SMichael Lange switch (s->nd) { 111f7320561SMichael Lange case 0: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed faces in Fluent files are not supported"); 112f7320561SMichael Lange case 2: numEntries = 2 + 2; break; /* linear */ 113f7320561SMichael Lange case 3: numEntries = 2 + 3; break; /* triangular */ 114f7320561SMichael Lange case 4: numEntries = 2 + 4; break; /* quadrilateral */ 115f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file"); 116f7320561SMichael Lange } 117f7320561SMichael Lange numFaces = s->last-s->first + 1; 118f7320561SMichael Lange ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr); 119f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 120f7320561SMichael Lange for (e = 0; e < numEntries; e++) { 121f7320561SMichael Lange ierr = PetscViewerRead(viewer, buffer, 1, PETSC_STRING);CHKERRQ(ierr); 122f7320561SMichael Lange snum = sscanf(buffer, "%x", &(((PetscInt*)s->data)[f*numEntries + e])); 123f7320561SMichael Lange if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file"); 124f7320561SMichael Lange } 125f7320561SMichael Lange } 126f7320561SMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1272f0bd6dcSMichael Lange } 1282f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr); 1292f0bd6dcSMichael Lange 1302f0bd6dcSMichael Lange } else { /* Unknown section type */ 1312f0bd6dcSMichael Lange PetscInt depth = 1; 1322f0bd6dcSMichael Lange do { 1332f0bd6dcSMichael Lange /* Match parentheses when parsing unknown sections */ 1342f0bd6dcSMichael Lange do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, PETSC_CHAR);CHKERRQ(ierr);} 1352f0bd6dcSMichael Lange while (buffer[0] != '(' && buffer[0] != ')'); 1362f0bd6dcSMichael Lange if (buffer[0] == '(') depth++; 1372f0bd6dcSMichael Lange if (buffer[0] == ')') depth--; 1382f0bd6dcSMichael Lange } while (depth > 0); 1392f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr); 1402f0bd6dcSMichael Lange } 1412f0bd6dcSMichael Lange PetscFunctionReturn(0); 1422f0bd6dcSMichael Lange } 1432f0bd6dcSMichael Lange 1442f0bd6dcSMichael Lange #undef __FUNCT__ 1452f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent" 1462f0bd6dcSMichael Lange /*@C 1472f0bd6dcSMichael Lange DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file. 1482f0bd6dcSMichael Lange 1492f0bd6dcSMichael Lange Collective on comm 1502f0bd6dcSMichael Lange 1512f0bd6dcSMichael Lange Input Parameters: 1522f0bd6dcSMichael Lange + comm - The MPI communicator 1532f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file 1542f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh 1552f0bd6dcSMichael Lange 1562f0bd6dcSMichael Lange Output Parameter: 1572f0bd6dcSMichael Lange . dm - The DM object representing the mesh 1582f0bd6dcSMichael Lange 1592f0bd6dcSMichael Lange Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm 1602f0bd6dcSMichael Lange 1612f0bd6dcSMichael Lange Level: beginner 1622f0bd6dcSMichael Lange 1632f0bd6dcSMichael Lange .keywords: mesh, fluent, case 1642f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate() 1652f0bd6dcSMichael Lange @*/ 1662f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm) 1672f0bd6dcSMichael Lange { 1682f0bd6dcSMichael Lange PetscMPIInt rank; 169f7320561SMichael Lange PetscInt c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE; 170f7320561SMichael Lange PetscInt numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE; 171f7320561SMichael Lange PetscInt *faces = NULL, *cellVertices; 172*3f6dc66eSMichael Lange PetscInt d, coordSize; 173*3f6dc66eSMichael Lange PetscScalar *coords, *coordsIn = NULL; 174*3f6dc66eSMichael Lange PetscSection coordSection; 175*3f6dc66eSMichael Lange Vec coordinates; 1762f0bd6dcSMichael Lange PetscErrorCode ierr; 1772f0bd6dcSMichael Lange 1782f0bd6dcSMichael Lange PetscFunctionBegin; 1792f0bd6dcSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1802f0bd6dcSMichael Lange 1812f0bd6dcSMichael Lange if (!rank) { 1822f0bd6dcSMichael Lange FluentSection s; 1832f0bd6dcSMichael Lange do { 1842f0bd6dcSMichael Lange ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr); 1852f0bd6dcSMichael Lange if (s.index == 2) { /* Dimension */ 186f7320561SMichael Lange dim = s.nd; 1872f0bd6dcSMichael Lange 1882f0bd6dcSMichael Lange } else if (s.index == 10) { /* Vertices */ 189f7320561SMichael Lange if (s.zoneID == 0) numVertices = s.last; 190*3f6dc66eSMichael Lange else { 191*3f6dc66eSMichael Lange if (coordsIn) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Currently no support for multiple coordinate sets in Fluent files"); 192*3f6dc66eSMichael Lange coordsIn = s.data; 193*3f6dc66eSMichael Lange } 1942f0bd6dcSMichael Lange 1952f0bd6dcSMichael Lange } else if (s.index == 12) { /* Cells */ 196f7320561SMichael Lange if (s.zoneID == 0) numCells = s.last; 197f7320561SMichael Lange else { 198f7320561SMichael Lange switch (s.nd) { 199f7320561SMichael Lange case 0: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed elements in Fluent files are not supported"); 200f7320561SMichael Lange case 1: numCellVertices = 3; break; /* triangular */ 201f7320561SMichael Lange case 2: numCellVertices = 4; break; /* tetrahedral */ 202f7320561SMichael Lange case 3: numCellVertices = 4; break; /* quadrilateral */ 203f7320561SMichael Lange case 4: numCellVertices = 8; break; /* hexahedral */ 204f7320561SMichael Lange case 5: numCellVertices = 5; break; /* pyramid */ 205f7320561SMichael Lange case 6: numCellVertices = 6; break; /* wedge */ 206f7320561SMichael Lange default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown cell element-type in Fluent file"); 207f7320561SMichael Lange } 208f7320561SMichael Lange } 2092f0bd6dcSMichael Lange 2102f0bd6dcSMichael Lange } else if (s.index == 13) { /* Facets */ 211f7320561SMichael Lange if (s.zoneID == 0) { /* Header section */ 212f7320561SMichael Lange numFaces = s.last - s.first + 1; 213f7320561SMichael Lange } else { /* Data section */ 214f7320561SMichael Lange if (s.nd == 0 || (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices)) { 215f7320561SMichael Lange SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are currently not supported"); 216f7320561SMichael Lange } 217f7320561SMichael Lange if (numFaces == PETSC_DETERMINE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file"); 218f7320561SMichael Lange if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd; 219f7320561SMichael Lange numFaceEntries = numFaceVertices + 2; 220f7320561SMichael Lange if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);} 221f7320561SMichael Lange ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr); 222f7320561SMichael Lange ierr = PetscFree(s.data);CHKERRQ(ierr); 223f7320561SMichael Lange } 2242f0bd6dcSMichael Lange } 2252f0bd6dcSMichael Lange } while (s.index >= 0); 2262f0bd6dcSMichael Lange } 227f7320561SMichael Lange ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 228f7320561SMichael Lange if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension"); 229f7320561SMichael Lange 230f7320561SMichael Lange /* Allocate cell-vertex mesh */ 231f7320561SMichael Lange ierr = DMCreate(comm, dm);CHKERRQ(ierr); 232f7320561SMichael Lange ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 233f7320561SMichael Lange ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 234f7320561SMichael Lange ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr); 235f7320561SMichael Lange if (!rank) { 236f7320561SMichael Lange if (numCellVertices == PETSC_DETERMINE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type not specified in Fluent file"); 237f7320561SMichael Lange for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);} 238f7320561SMichael Lange } 239f7320561SMichael Lange ierr = DMSetUp(*dm);CHKERRQ(ierr); 240f7320561SMichael Lange 241f7320561SMichael Lange if (!rank) { 242f7320561SMichael Lange /* Derive cell-vertex list from face-vertex and face-cell maps */ 243f7320561SMichael Lange if (numCells == PETSC_DETERMINE || numCellVertices == PETSC_DETERMINE) { 244f7320561SMichael Lange SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Insufficent cell header information in Fluent file"); 245f7320561SMichael Lange } 246f7320561SMichael Lange ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr); 247f7320561SMichael Lange for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1; 248f7320561SMichael Lange for (f = 0; f < numFaces; f++) { 249f7320561SMichael Lange PetscInt *cell; 250f7320561SMichael Lange const PetscInt cl = faces[f*numFaceEntries + numFaceVertices]; 251f7320561SMichael Lange const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1]; 252f7320561SMichael Lange const PetscInt *face = &(faces[f*numFaceEntries]); 253f7320561SMichael Lange 254f7320561SMichael Lange if (cl > 0) { 255f7320561SMichael Lange cell = &(cellVertices[(cl-1) * numCellVertices]); 256f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 257f7320561SMichael Lange PetscBool found = PETSC_FALSE; 258f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 259f7320561SMichael Lange if (cell[c] < 0) break; 260f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 261f7320561SMichael Lange } 262f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 263f7320561SMichael Lange } 264f7320561SMichael Lange } 265f7320561SMichael Lange if (cr > 0) { 266f7320561SMichael Lange cell = &(cellVertices[(cr-1) * numCellVertices]); 267f7320561SMichael Lange for (v = 0; v < numFaceVertices; v++) { 268f7320561SMichael Lange PetscBool found = PETSC_FALSE; 269f7320561SMichael Lange for (c = 0; c < numCellVertices; c++) { 270f7320561SMichael Lange if (cell[c] < 0) break; 271f7320561SMichael Lange if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;} 272f7320561SMichael Lange } 273f7320561SMichael Lange if (!found) cell[c] = face[v]-1 + numCells; 274f7320561SMichael Lange } 275f7320561SMichael Lange } 276f7320561SMichael Lange } 277f7320561SMichael Lange for (c = 0; c < numCells; c++) { 278f7320561SMichael Lange ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr); 279f7320561SMichael Lange } 280f7320561SMichael Lange } 281f7320561SMichael Lange ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 282f7320561SMichael Lange ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 283f7320561SMichael Lange if (interpolate) { 284f7320561SMichael Lange DM idm = NULL; 285f7320561SMichael Lange 286f7320561SMichael Lange ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 287f7320561SMichael Lange ierr = DMDestroy(dm);CHKERRQ(ierr); 288f7320561SMichael Lange *dm = idm; 289f7320561SMichael Lange } 290f7320561SMichael Lange 291*3f6dc66eSMichael Lange /* Read coordinates */ 292*3f6dc66eSMichael Lange ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 293*3f6dc66eSMichael Lange ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 294*3f6dc66eSMichael Lange ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 295*3f6dc66eSMichael Lange ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 296*3f6dc66eSMichael Lange for (v = numCells; v < numCells+numVertices; ++v) { 297*3f6dc66eSMichael Lange ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 298*3f6dc66eSMichael Lange ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 299*3f6dc66eSMichael Lange } 300*3f6dc66eSMichael Lange ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 301*3f6dc66eSMichael Lange ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 302*3f6dc66eSMichael Lange ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 303*3f6dc66eSMichael Lange ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 304*3f6dc66eSMichael Lange ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 305*3f6dc66eSMichael Lange ierr = VecSetType(coordinates, VECSTANDARD);CHKERRQ(ierr); 306*3f6dc66eSMichael Lange ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 307*3f6dc66eSMichael Lange if (!rank) { 308*3f6dc66eSMichael Lange for (v = 0; v < numVertices; ++v) { 309*3f6dc66eSMichael Lange for (d = 0; d < dim; ++d) { 310*3f6dc66eSMichael Lange coords[v*dim+d] = coordsIn[v*dim+d]; 311*3f6dc66eSMichael Lange } 312*3f6dc66eSMichael Lange } 313*3f6dc66eSMichael Lange } 314*3f6dc66eSMichael Lange ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 315*3f6dc66eSMichael Lange ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 316*3f6dc66eSMichael Lange ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 317f7320561SMichael Lange if (!rank) { 318f7320561SMichael Lange ierr = PetscFree(cellVertices);CHKERRQ(ierr); 319f7320561SMichael Lange ierr = PetscFree(faces);CHKERRQ(ierr); 320*3f6dc66eSMichael Lange ierr = PetscFree(coordsIn);CHKERRQ(ierr); 321f7320561SMichael Lange } 3222f0bd6dcSMichael Lange PetscFunctionReturn(0); 3232f0bd6dcSMichael Lange } 324