xref: /petsc/src/dm/impls/plex/plexfluent.c (revision f7320561e3b3aa527dec50a409aa9784d53fbd62)
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");
782f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
792f0bd6dcSMichael Lange 
802f0bd6dcSMichael Lange   } else if (s->index == 12) {   /* Cells */
812f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
822f0bd6dcSMichael Lange     snum = sscanf(buffer, "(%x", &(s->zoneID));
832f0bd6dcSMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
842f0bd6dcSMichael Lange     if (s->zoneID == 0) {  /* Header section */
852f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
862f0bd6dcSMichael Lange       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
87*f7320561SMichael Lange     } else {               /* Data section */
882f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
892f0bd6dcSMichael Lange       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
902f0bd6dcSMichael Lange     }
912f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
922f0bd6dcSMichael Lange 
932f0bd6dcSMichael Lange   } else if (s->index == 13) {   /* Faces */
942f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
952f0bd6dcSMichael Lange     snum = sscanf(buffer, "(%x", &(s->zoneID));
962f0bd6dcSMichael Lange     if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
972f0bd6dcSMichael Lange     if (s->zoneID == 0) {  /* Header section */
982f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d)", &(s->zoneID), &(s->first), &(s->last), &(s->nd));
992f0bd6dcSMichael Lange       if (snum != 4) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
100*f7320561SMichael Lange     } else {               /* Data section */
101*f7320561SMichael Lange       PetscInt f, e, numEntries, numFaces;
1022f0bd6dcSMichael Lange       snum = sscanf(buffer, "(%x %x %x %d %d)", &(s->zoneID), &(s->first), &(s->last), &(s->type), &(s->nd));
1032f0bd6dcSMichael Lange       if (snum != 5) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
104*f7320561SMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '(');CHKERRQ(ierr);
105*f7320561SMichael Lange       switch (s->nd) {
106*f7320561SMichael Lange       case 0: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed faces in Fluent files are not supported");
107*f7320561SMichael Lange       case 2: numEntries = 2 + 2; break;  /* linear */
108*f7320561SMichael Lange       case 3: numEntries = 2 + 3; break;  /* triangular */
109*f7320561SMichael Lange       case 4: numEntries = 2 + 4; break;  /* quadrilateral */
110*f7320561SMichael Lange       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown face type in Fluent file");
111*f7320561SMichael Lange       }
112*f7320561SMichael Lange       numFaces = s->last-s->first + 1;
113*f7320561SMichael Lange       ierr = PetscMalloc1(numEntries*numFaces, (PetscInt**)&s->data);CHKERRQ(ierr);
114*f7320561SMichael Lange       for (f = 0; f < numFaces; f++) {
115*f7320561SMichael Lange         for (e = 0; e < numEntries; e++) {
116*f7320561SMichael Lange           ierr = PetscViewerRead(viewer, buffer, 1, PETSC_STRING);CHKERRQ(ierr);
117*f7320561SMichael Lange           snum = sscanf(buffer, "%x", &(((PetscInt*)s->data)[f*numEntries + e]));
118*f7320561SMichael Lange           if (snum != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "File is not a valid Fluent file");
119*f7320561SMichael Lange         }
120*f7320561SMichael Lange       }
121*f7320561SMichael Lange       ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1222f0bd6dcSMichael Lange     }
1232f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, ')');CHKERRQ(ierr);
1242f0bd6dcSMichael Lange 
1252f0bd6dcSMichael Lange   } else {                       /* Unknown section type */
1262f0bd6dcSMichael Lange     PetscInt depth = 1;
1272f0bd6dcSMichael Lange     do {
1282f0bd6dcSMichael Lange       /* Match parentheses when parsing unknown sections */
1292f0bd6dcSMichael Lange       do {ierr = PetscViewerRead(viewer, &(buffer[0]), 1, PETSC_CHAR);CHKERRQ(ierr);}
1302f0bd6dcSMichael Lange       while (buffer[0] != '(' && buffer[0] != ')');
1312f0bd6dcSMichael Lange       if (buffer[0] == '(') depth++;
1322f0bd6dcSMichael Lange       if (buffer[0] == ')') depth--;
1332f0bd6dcSMichael Lange     } while (depth > 0);
1342f0bd6dcSMichael Lange     ierr = DMPlexCreateFluent_ReadString(viewer, buffer, '\n');CHKERRQ(ierr);
1352f0bd6dcSMichael Lange   }
1362f0bd6dcSMichael Lange   PetscFunctionReturn(0);
1372f0bd6dcSMichael Lange }
1382f0bd6dcSMichael Lange 
1392f0bd6dcSMichael Lange #undef __FUNCT__
1402f0bd6dcSMichael Lange #define __FUNCT__ "DMPlexCreateFluent"
1412f0bd6dcSMichael Lange /*@C
1422f0bd6dcSMichael Lange   DMPlexCreateFluent - Create a DMPlex mesh from a Fluent mesh file.
1432f0bd6dcSMichael Lange 
1442f0bd6dcSMichael Lange   Collective on comm
1452f0bd6dcSMichael Lange 
1462f0bd6dcSMichael Lange   Input Parameters:
1472f0bd6dcSMichael Lange + comm  - The MPI communicator
1482f0bd6dcSMichael Lange . viewer - The Viewer associated with a Fluent mesh file
1492f0bd6dcSMichael Lange - interpolate - Create faces and edges in the mesh
1502f0bd6dcSMichael Lange 
1512f0bd6dcSMichael Lange   Output Parameter:
1522f0bd6dcSMichael Lange . dm  - The DM object representing the mesh
1532f0bd6dcSMichael Lange 
1542f0bd6dcSMichael Lange   Note: http://aerojet.engr.ucdavis.edu/fluenthelp/html/ug/node1490.htm
1552f0bd6dcSMichael Lange 
1562f0bd6dcSMichael Lange   Level: beginner
1572f0bd6dcSMichael Lange 
1582f0bd6dcSMichael Lange .keywords: mesh, fluent, case
1592f0bd6dcSMichael Lange .seealso: DMPLEX, DMCreate()
1602f0bd6dcSMichael Lange @*/
1612f0bd6dcSMichael Lange PetscErrorCode DMPlexCreateFluent(MPI_Comm comm, PetscViewer viewer, PetscBool interpolate, DM *dm)
1622f0bd6dcSMichael Lange {
1632f0bd6dcSMichael Lange   PetscMPIInt    rank;
164*f7320561SMichael Lange   PetscInt       c, f, v, dim = PETSC_DETERMINE, numCells = 0, numVertices = 0, numCellVertices = PETSC_DETERMINE;
165*f7320561SMichael Lange   PetscInt       numFaces = PETSC_DETERMINE, numFaceEntries = PETSC_DETERMINE, numFaceVertices = PETSC_DETERMINE;
166*f7320561SMichael Lange   PetscInt      *faces = NULL, *cellVertices;
1672f0bd6dcSMichael Lange   PetscErrorCode ierr;
1682f0bd6dcSMichael Lange 
1692f0bd6dcSMichael Lange   PetscFunctionBegin;
1702f0bd6dcSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1712f0bd6dcSMichael Lange 
1722f0bd6dcSMichael Lange   if (!rank) {
1732f0bd6dcSMichael Lange     FluentSection s;
1742f0bd6dcSMichael Lange     do {
1752f0bd6dcSMichael Lange       ierr = DMPlexCreateFluent_ReadSection(viewer, &s);CHKERRQ(ierr);
1762f0bd6dcSMichael Lange       if (s.index == 2) {                 /* Dimension */
177*f7320561SMichael Lange         dim = s.nd;
1782f0bd6dcSMichael Lange 
1792f0bd6dcSMichael Lange       } else if (s.index == 10) {         /* Vertices */
180*f7320561SMichael Lange         if (s.zoneID == 0) numVertices = s.last;
1812f0bd6dcSMichael Lange 
1822f0bd6dcSMichael Lange       } else if (s.index == 12) {         /* Cells */
183*f7320561SMichael Lange         if (s.zoneID == 0) numCells = s.last;
184*f7320561SMichael Lange         else {
185*f7320561SMichael Lange           switch (s.nd) {
186*f7320561SMichael Lange           case 0: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed elements in Fluent files are not supported");
187*f7320561SMichael Lange           case 1: numCellVertices = 3; break;  /* triangular */
188*f7320561SMichael Lange           case 2: numCellVertices = 4; break;  /* tetrahedral */
189*f7320561SMichael Lange           case 3: numCellVertices = 4; break;  /* quadrilateral */
190*f7320561SMichael Lange           case 4: numCellVertices = 8; break;  /* hexahedral */
191*f7320561SMichael Lange           case 5: numCellVertices = 5; break;  /* pyramid */
192*f7320561SMichael Lange           case 6: numCellVertices = 6; break;  /* wedge */
193*f7320561SMichael Lange           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unknown cell element-type in Fluent file");
194*f7320561SMichael Lange           }
195*f7320561SMichael Lange         }
1962f0bd6dcSMichael Lange 
1972f0bd6dcSMichael Lange       } else if (s.index == 13) {         /* Facets */
198*f7320561SMichael Lange         if (s.zoneID == 0) {  /* Header section */
199*f7320561SMichael Lange           numFaces = s.last - s.first + 1;
200*f7320561SMichael Lange         } else {              /* Data section */
201*f7320561SMichael Lange           if (s.nd == 0 || (numFaceVertices != PETSC_DETERMINE && s.nd != numFaceVertices)) {
202*f7320561SMichael Lange             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Mixed facets in Fluent files are currently not supported");
203*f7320561SMichael Lange           }
204*f7320561SMichael Lange           if (numFaces == PETSC_DETERMINE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No header section for facets in Fluent file");
205*f7320561SMichael Lange           if (numFaceVertices == PETSC_DETERMINE) numFaceVertices = s.nd;
206*f7320561SMichael Lange           numFaceEntries = numFaceVertices + 2;
207*f7320561SMichael Lange           if (!faces) {ierr = PetscMalloc1(numFaces*numFaceEntries, &faces);CHKERRQ(ierr);}
208*f7320561SMichael Lange           ierr = PetscMemcpy(&(faces[(s.first-1)*numFaceEntries]), s.data, (s.last-s.first+1)*numFaceEntries*sizeof(PetscInt));CHKERRQ(ierr);
209*f7320561SMichael Lange           ierr = PetscFree(s.data);CHKERRQ(ierr);
210*f7320561SMichael Lange         }
2112f0bd6dcSMichael Lange       }
2122f0bd6dcSMichael Lange     } while (s.index >= 0);
2132f0bd6dcSMichael Lange   }
214*f7320561SMichael Lange   ierr = MPI_Bcast(&dim, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
215*f7320561SMichael Lange   if (dim < 0) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Fluent file does not include dimension");
216*f7320561SMichael Lange 
217*f7320561SMichael Lange   /* Allocate cell-vertex mesh */
218*f7320561SMichael Lange   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
219*f7320561SMichael Lange   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
220*f7320561SMichael Lange   ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr);
221*f7320561SMichael Lange   ierr = DMPlexSetChart(*dm, 0, numCells + numVertices);CHKERRQ(ierr);
222*f7320561SMichael Lange   if (!rank) {
223*f7320561SMichael Lange     if (numCellVertices == PETSC_DETERMINE) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Cell type not specified in Fluent file");
224*f7320561SMichael Lange     for (c = 0; c < numCells; ++c) {ierr = DMPlexSetConeSize(*dm, c, numCellVertices);CHKERRQ(ierr);}
225*f7320561SMichael Lange   }
226*f7320561SMichael Lange   ierr = DMSetUp(*dm);CHKERRQ(ierr);
227*f7320561SMichael Lange 
228*f7320561SMichael Lange   if (!rank) {
229*f7320561SMichael Lange     /* Derive cell-vertex list from face-vertex and face-cell maps */
230*f7320561SMichael Lange     if (numCells == PETSC_DETERMINE || numCellVertices == PETSC_DETERMINE) {
231*f7320561SMichael Lange       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Insufficent cell header information in Fluent file");
232*f7320561SMichael Lange     }
233*f7320561SMichael Lange     ierr = PetscMalloc1(numCells*numCellVertices, &cellVertices);CHKERRQ(ierr);
234*f7320561SMichael Lange     for (c = 0; c < numCells*numCellVertices; c++) cellVertices[c] = -1;
235*f7320561SMichael Lange     for (f = 0; f < numFaces; f++) {
236*f7320561SMichael Lange       PetscInt *cell;
237*f7320561SMichael Lange       const PetscInt cl = faces[f*numFaceEntries + numFaceVertices];
238*f7320561SMichael Lange       const PetscInt cr = faces[f*numFaceEntries + numFaceVertices + 1];
239*f7320561SMichael Lange       const PetscInt *face = &(faces[f*numFaceEntries]);
240*f7320561SMichael Lange 
241*f7320561SMichael Lange       if (cl > 0) {
242*f7320561SMichael Lange         cell = &(cellVertices[(cl-1) * numCellVertices]);
243*f7320561SMichael Lange         for (v = 0; v < numFaceVertices; v++) {
244*f7320561SMichael Lange           PetscBool found = PETSC_FALSE;
245*f7320561SMichael Lange           for (c = 0; c < numCellVertices; c++) {
246*f7320561SMichael Lange             if (cell[c] < 0) break;
247*f7320561SMichael Lange             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
248*f7320561SMichael Lange           }
249*f7320561SMichael Lange           if (!found) cell[c] = face[v]-1 + numCells;
250*f7320561SMichael Lange         }
251*f7320561SMichael Lange       }
252*f7320561SMichael Lange       if (cr > 0) {
253*f7320561SMichael Lange         cell = &(cellVertices[(cr-1) * numCellVertices]);
254*f7320561SMichael Lange         for (v = 0; v < numFaceVertices; v++) {
255*f7320561SMichael Lange           PetscBool found = PETSC_FALSE;
256*f7320561SMichael Lange           for (c = 0; c < numCellVertices; c++) {
257*f7320561SMichael Lange             if (cell[c] < 0) break;
258*f7320561SMichael Lange             if (cell[c] == face[v]-1 + numCells) {found = PETSC_TRUE; break;}
259*f7320561SMichael Lange           }
260*f7320561SMichael Lange           if (!found) cell[c] = face[v]-1 + numCells;
261*f7320561SMichael Lange         }
262*f7320561SMichael Lange       }
263*f7320561SMichael Lange     }
264*f7320561SMichael Lange     for (c = 0; c < numCells; c++) {
265*f7320561SMichael Lange       ierr = DMPlexSetCone(*dm, c, &(cellVertices[c*numCellVertices]));CHKERRQ(ierr);
266*f7320561SMichael Lange     }
267*f7320561SMichael Lange   }
268*f7320561SMichael Lange   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
269*f7320561SMichael Lange   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
270*f7320561SMichael Lange   if (interpolate) {
271*f7320561SMichael Lange     DM idm = NULL;
272*f7320561SMichael Lange 
273*f7320561SMichael Lange     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
274*f7320561SMichael Lange     ierr = DMDestroy(dm);CHKERRQ(ierr);
275*f7320561SMichael Lange     *dm  = idm;
276*f7320561SMichael Lange   }
277*f7320561SMichael Lange 
278*f7320561SMichael Lange   if (!rank) {
279*f7320561SMichael Lange     ierr = PetscFree(cellVertices);CHKERRQ(ierr);
280*f7320561SMichael Lange     ierr = PetscFree(faces);CHKERRQ(ierr);
281*f7320561SMichael Lange   }
2822f0bd6dcSMichael Lange   PetscFunctionReturn(0);
2832f0bd6dcSMichael Lange }
284