xref: /petsc/src/dm/impls/plex/cgns/plexcgns2.c (revision 1fca310d5bbec1c360cd9da8976f4e0178ab9137)
15f34f2dcSJed Brown #define PETSCDM_DLL
25f34f2dcSJed Brown #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
35f34f2dcSJed Brown #include <petsc/private/viewercgnsimpl.h>
45f34f2dcSJed Brown 
55f34f2dcSJed Brown #include <pcgnslib.h>
65f34f2dcSJed Brown #include <cgns_io.h>
75f34f2dcSJed Brown 
85f34f2dcSJed Brown #if !defined(CGNS_ENUMT)
95f34f2dcSJed Brown   #define CGNS_ENUMT(a) a
105f34f2dcSJed Brown #endif
115f34f2dcSJed Brown #if !defined(CGNS_ENUMV)
125f34f2dcSJed Brown   #define CGNS_ENUMV(a) a
135f34f2dcSJed Brown #endif
145f34f2dcSJed Brown 
15d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
16d71ae5a4SJacob Faibussowitsch {
175f34f2dcSJed Brown   PetscMPIInt rank;
185f34f2dcSJed Brown   int         cgid = -1;
195f34f2dcSJed Brown 
205f34f2dcSJed Brown   PetscFunctionBegin;
214f572ea9SToby Isaac   PetscAssertPointer(filename, 2);
225f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
235f34f2dcSJed Brown   if (rank == 0) {
245f34f2dcSJed Brown     PetscCallCGNS(cg_open(filename, CG_MODE_READ, &cgid));
255f34f2dcSJed Brown     PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename);
265f34f2dcSJed Brown   }
275f34f2dcSJed Brown   PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm));
285f34f2dcSJed Brown   if (rank == 0) PetscCallCGNS(cg_close(cgid));
293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305f34f2dcSJed Brown }
315f34f2dcSJed Brown 
32d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
33d71ae5a4SJacob Faibussowitsch {
345f34f2dcSJed Brown   PetscMPIInt  num_proc, rank;
355f34f2dcSJed Brown   DM           cdm;
365f34f2dcSJed Brown   DMLabel      label;
375f34f2dcSJed Brown   PetscSection coordSection;
385f34f2dcSJed Brown   Vec          coordinates;
395f34f2dcSJed Brown   PetscScalar *coords;
405f34f2dcSJed Brown   PetscInt    *cellStart, *vertStart, v;
415f34f2dcSJed Brown   PetscInt     labelIdRange[2], labelId;
425f34f2dcSJed Brown   /* Read from file */
435f34f2dcSJed Brown   char basename[CGIO_MAX_NAME_LENGTH + 1];
445f34f2dcSJed Brown   char buffer[CGIO_MAX_NAME_LENGTH + 1];
455f34f2dcSJed Brown   int  dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0;
465f34f2dcSJed Brown   int  nzones = 0;
475f34f2dcSJed Brown 
485f34f2dcSJed Brown   PetscFunctionBegin;
495f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_rank(comm, &rank));
505f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
515f34f2dcSJed Brown   PetscCall(DMCreate(comm, dm));
525f34f2dcSJed Brown   PetscCall(DMSetType(*dm, DMPLEX));
535f34f2dcSJed Brown 
545f34f2dcSJed Brown   /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */
555f34f2dcSJed Brown   if (rank == 0) {
565f34f2dcSJed Brown     int nbases, z;
575f34f2dcSJed Brown 
585f34f2dcSJed Brown     PetscCallCGNS(cg_nbases(cgid, &nbases));
595f34f2dcSJed Brown     PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases);
605f34f2dcSJed Brown     PetscCallCGNS(cg_base_read(cgid, 1, basename, &dim, &physDim));
615f34f2dcSJed Brown     PetscCallCGNS(cg_nzones(cgid, 1, &nzones));
625f34f2dcSJed Brown     PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart));
635f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
645f34f2dcSJed Brown       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
655f34f2dcSJed Brown 
665f34f2dcSJed Brown       PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
675f34f2dcSJed Brown       numVertices += sizes[0];
685f34f2dcSJed Brown       numCells += sizes[1];
695f34f2dcSJed Brown       cellStart[z] += sizes[1] + cellStart[z - 1];
705f34f2dcSJed Brown       vertStart[z] += sizes[0] + vertStart[z - 1];
715f34f2dcSJed Brown     }
72ad540459SPierre Jolivet     for (z = 1; z <= nzones; ++z) vertStart[z] += numCells;
735f34f2dcSJed Brown     coordDim = dim;
745f34f2dcSJed Brown   }
755f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm));
765f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
775f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm));
785f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm));
795f34f2dcSJed Brown 
805f34f2dcSJed Brown   PetscCall(PetscObjectSetName((PetscObject)*dm, basename));
815f34f2dcSJed Brown   PetscCall(DMSetDimension(*dm, dim));
825f34f2dcSJed Brown   PetscCall(DMCreateLabel(*dm, "celltype"));
835f34f2dcSJed Brown   PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
845f34f2dcSJed Brown 
855f34f2dcSJed Brown   /* Read zone information */
865f34f2dcSJed Brown   if (rank == 0) {
875f34f2dcSJed Brown     int z, c, c_loc;
885f34f2dcSJed Brown 
895f34f2dcSJed Brown     /* Read the cell set connectivity table and build mesh topology
905f34f2dcSJed Brown        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
915f34f2dcSJed Brown     /* First set sizes */
925f34f2dcSJed Brown     for (z = 1, c = 0; z <= nzones; ++z) {
935f34f2dcSJed Brown       CGNS_ENUMT(ZoneType_t) zonetype;
945f34f2dcSJed Brown       int nsections;
955f34f2dcSJed Brown       CGNS_ENUMT(ElementType_t) cellType;
965f34f2dcSJed Brown       cgsize_t       start, end;
975f34f2dcSJed Brown       int            nbndry, parentFlag;
985f34f2dcSJed Brown       PetscInt       numCorners;
995f34f2dcSJed Brown       DMPolytopeType ctype;
1005f34f2dcSJed Brown 
1015f34f2dcSJed Brown       PetscCallCGNS(cg_zone_type(cgid, 1, z, &zonetype));
1025f34f2dcSJed Brown       PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS");
1035f34f2dcSJed Brown       PetscCallCGNS(cg_nsections(cgid, 1, z, &nsections));
1045f34f2dcSJed Brown       PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections);
1055f34f2dcSJed Brown       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
1065f34f2dcSJed Brown       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
1075f34f2dcSJed Brown       if (cellType == CGNS_ENUMV(MIXED)) {
1085f34f2dcSJed Brown         cgsize_t elementDataSize, *elements;
1095f34f2dcSJed Brown         PetscInt off;
1105f34f2dcSJed Brown 
1115f34f2dcSJed Brown         PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
1125f34f2dcSJed Brown         PetscCall(PetscMalloc1(elementDataSize, &elements));
1135f34f2dcSJed Brown         PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
1145f34f2dcSJed Brown         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
1155f34f2dcSJed Brown           switch (elements[off]) {
1169371c9d4SSatish Balay           case CGNS_ENUMV(BAR_2):
1179371c9d4SSatish Balay             numCorners = 2;
1189371c9d4SSatish Balay             ctype      = DM_POLYTOPE_SEGMENT;
1199371c9d4SSatish Balay             break;
1209371c9d4SSatish Balay           case CGNS_ENUMV(TRI_3):
1219371c9d4SSatish Balay             numCorners = 3;
1229371c9d4SSatish Balay             ctype      = DM_POLYTOPE_TRIANGLE;
1239371c9d4SSatish Balay             break;
1249371c9d4SSatish Balay           case CGNS_ENUMV(QUAD_4):
1259371c9d4SSatish Balay             numCorners = 4;
1269371c9d4SSatish Balay             ctype      = DM_POLYTOPE_QUADRILATERAL;
1279371c9d4SSatish Balay             break;
1289371c9d4SSatish Balay           case CGNS_ENUMV(TETRA_4):
1299371c9d4SSatish Balay             numCorners = 4;
1309371c9d4SSatish Balay             ctype      = DM_POLYTOPE_TETRAHEDRON;
1319371c9d4SSatish Balay             break;
1329371c9d4SSatish Balay           case CGNS_ENUMV(PENTA_6):
1339371c9d4SSatish Balay             numCorners = 6;
1349371c9d4SSatish Balay             ctype      = DM_POLYTOPE_TRI_PRISM;
1359371c9d4SSatish Balay             break;
1369371c9d4SSatish Balay           case CGNS_ENUMV(HEXA_8):
1379371c9d4SSatish Balay             numCorners = 8;
1389371c9d4SSatish Balay             ctype      = DM_POLYTOPE_HEXAHEDRON;
1399371c9d4SSatish Balay             break;
140d71ae5a4SJacob Faibussowitsch           default:
141d71ae5a4SJacob Faibussowitsch             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)elements[off]);
1425f34f2dcSJed Brown           }
1435f34f2dcSJed Brown           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
1445f34f2dcSJed Brown           PetscCall(DMPlexSetCellType(*dm, c, ctype));
1455f34f2dcSJed Brown           off += numCorners + 1;
1465f34f2dcSJed Brown         }
1475f34f2dcSJed Brown         PetscCall(PetscFree(elements));
1485f34f2dcSJed Brown       } else {
1495f34f2dcSJed Brown         switch (cellType) {
1509371c9d4SSatish Balay         case CGNS_ENUMV(BAR_2):
1519371c9d4SSatish Balay           numCorners = 2;
1529371c9d4SSatish Balay           ctype      = DM_POLYTOPE_SEGMENT;
1539371c9d4SSatish Balay           break;
1549371c9d4SSatish Balay         case CGNS_ENUMV(TRI_3):
1559371c9d4SSatish Balay           numCorners = 3;
1569371c9d4SSatish Balay           ctype      = DM_POLYTOPE_TRIANGLE;
1579371c9d4SSatish Balay           break;
1589371c9d4SSatish Balay         case CGNS_ENUMV(QUAD_4):
1599371c9d4SSatish Balay           numCorners = 4;
1609371c9d4SSatish Balay           ctype      = DM_POLYTOPE_QUADRILATERAL;
1619371c9d4SSatish Balay           break;
1629371c9d4SSatish Balay         case CGNS_ENUMV(TETRA_4):
1639371c9d4SSatish Balay           numCorners = 4;
1649371c9d4SSatish Balay           ctype      = DM_POLYTOPE_TETRAHEDRON;
1659371c9d4SSatish Balay           break;
1669371c9d4SSatish Balay         case CGNS_ENUMV(PENTA_6):
1679371c9d4SSatish Balay           numCorners = 6;
1689371c9d4SSatish Balay           ctype      = DM_POLYTOPE_TRI_PRISM;
1699371c9d4SSatish Balay           break;
1709371c9d4SSatish Balay         case CGNS_ENUMV(HEXA_8):
1719371c9d4SSatish Balay           numCorners = 8;
1729371c9d4SSatish Balay           ctype      = DM_POLYTOPE_HEXAHEDRON;
1739371c9d4SSatish Balay           break;
174d71ae5a4SJacob Faibussowitsch         default:
175d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)cellType);
1765f34f2dcSJed Brown         }
1775f34f2dcSJed Brown         for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
1785f34f2dcSJed Brown           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
1795f34f2dcSJed Brown           PetscCall(DMPlexSetCellType(*dm, c, ctype));
1805f34f2dcSJed Brown         }
1815f34f2dcSJed Brown       }
1825f34f2dcSJed Brown     }
18348a46eb9SPierre Jolivet     for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
1845f34f2dcSJed Brown   }
1855f34f2dcSJed Brown 
1865f34f2dcSJed Brown   PetscCall(DMSetUp(*dm));
1875f34f2dcSJed Brown 
1885f34f2dcSJed Brown   PetscCall(DMCreateLabel(*dm, "zone"));
1895f34f2dcSJed Brown   if (rank == 0) {
1905f34f2dcSJed Brown     int z, c, c_loc, v_loc;
1915f34f2dcSJed Brown 
1925f34f2dcSJed Brown     PetscCall(DMGetLabel(*dm, "zone", &label));
1935f34f2dcSJed Brown     for (z = 1, c = 0; z <= nzones; ++z) {
1945f34f2dcSJed Brown       CGNS_ENUMT(ElementType_t) cellType;
1955f34f2dcSJed Brown       cgsize_t  elementDataSize, *elements, start, end;
1965f34f2dcSJed Brown       int       nbndry, parentFlag;
1975f34f2dcSJed Brown       PetscInt *cone, numc, numCorners, maxCorners = 27;
1985f34f2dcSJed Brown 
1995f34f2dcSJed Brown       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
2005f34f2dcSJed Brown       numc = end - start;
2015f34f2dcSJed Brown       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
2025f34f2dcSJed Brown       PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
2035f34f2dcSJed Brown       PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone));
2045f34f2dcSJed Brown       PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
2055f34f2dcSJed Brown       if (cellType == CGNS_ENUMV(MIXED)) {
2065f34f2dcSJed Brown         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
2075f34f2dcSJed Brown         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
2085f34f2dcSJed Brown           switch (elements[v]) {
209d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(BAR_2):
210d71ae5a4SJacob Faibussowitsch             numCorners = 2;
211d71ae5a4SJacob Faibussowitsch             break;
212d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(TRI_3):
213d71ae5a4SJacob Faibussowitsch             numCorners = 3;
214d71ae5a4SJacob Faibussowitsch             break;
215d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(QUAD_4):
216d71ae5a4SJacob Faibussowitsch             numCorners = 4;
217d71ae5a4SJacob Faibussowitsch             break;
218d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(TETRA_4):
219d71ae5a4SJacob Faibussowitsch             numCorners = 4;
220d71ae5a4SJacob Faibussowitsch             break;
221d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(PENTA_6):
222d71ae5a4SJacob Faibussowitsch             numCorners = 6;
223d71ae5a4SJacob Faibussowitsch             break;
224d71ae5a4SJacob Faibussowitsch           case CGNS_ENUMV(HEXA_8):
225d71ae5a4SJacob Faibussowitsch             numCorners = 8;
226d71ae5a4SJacob Faibussowitsch             break;
227d71ae5a4SJacob Faibussowitsch           default:
228d71ae5a4SJacob Faibussowitsch             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)elements[v]);
2295f34f2dcSJed Brown           }
2305f34f2dcSJed Brown           ++v;
231ad540459SPierre Jolivet           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
2325f34f2dcSJed Brown           PetscCall(DMPlexReorderCell(*dm, c, cone));
2335f34f2dcSJed Brown           PetscCall(DMPlexSetCone(*dm, c, cone));
2345f34f2dcSJed Brown           PetscCall(DMLabelSetValue(label, c, z));
2355f34f2dcSJed Brown         }
2365f34f2dcSJed Brown       } else {
2375f34f2dcSJed Brown         switch (cellType) {
238d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(BAR_2):
239d71ae5a4SJacob Faibussowitsch           numCorners = 2;
240d71ae5a4SJacob Faibussowitsch           break;
241d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(TRI_3):
242d71ae5a4SJacob Faibussowitsch           numCorners = 3;
243d71ae5a4SJacob Faibussowitsch           break;
244d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(QUAD_4):
245d71ae5a4SJacob Faibussowitsch           numCorners = 4;
246d71ae5a4SJacob Faibussowitsch           break;
247d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(TETRA_4):
248d71ae5a4SJacob Faibussowitsch           numCorners = 4;
249d71ae5a4SJacob Faibussowitsch           break;
250d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(PENTA_6):
251d71ae5a4SJacob Faibussowitsch           numCorners = 6;
252d71ae5a4SJacob Faibussowitsch           break;
253d71ae5a4SJacob Faibussowitsch         case CGNS_ENUMV(HEXA_8):
254d71ae5a4SJacob Faibussowitsch           numCorners = 8;
255d71ae5a4SJacob Faibussowitsch           break;
256d71ae5a4SJacob Faibussowitsch         default:
257d71ae5a4SJacob Faibussowitsch           SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int)cellType);
2585f34f2dcSJed Brown         }
2595f34f2dcSJed Brown         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
2605f34f2dcSJed Brown         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
261ad540459SPierre Jolivet           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1;
2625f34f2dcSJed Brown           PetscCall(DMPlexReorderCell(*dm, c, cone));
2635f34f2dcSJed Brown           PetscCall(DMPlexSetCone(*dm, c, cone));
2645f34f2dcSJed Brown           PetscCall(DMLabelSetValue(label, c, z));
2655f34f2dcSJed Brown         }
2665f34f2dcSJed Brown       }
2675f34f2dcSJed Brown       PetscCall(PetscFree2(elements, cone));
2685f34f2dcSJed Brown     }
2695f34f2dcSJed Brown   }
2705f34f2dcSJed Brown 
2715f34f2dcSJed Brown   PetscCall(DMPlexSymmetrize(*dm));
2725f34f2dcSJed Brown   PetscCall(DMPlexStratify(*dm));
2735f34f2dcSJed Brown   if (interpolate) {
2745f34f2dcSJed Brown     DM idm;
2755f34f2dcSJed Brown 
2765f34f2dcSJed Brown     PetscCall(DMPlexInterpolate(*dm, &idm));
2775f34f2dcSJed Brown     PetscCall(DMDestroy(dm));
2785f34f2dcSJed Brown     *dm = idm;
2795f34f2dcSJed Brown   }
2805f34f2dcSJed Brown 
2815f34f2dcSJed Brown   /* Read coordinates */
2825f34f2dcSJed Brown   PetscCall(DMSetCoordinateDim(*dm, coordDim));
2835f34f2dcSJed Brown   PetscCall(DMGetCoordinateDM(*dm, &cdm));
2845f34f2dcSJed Brown   PetscCall(DMGetLocalSection(cdm, &coordSection));
2855f34f2dcSJed Brown   PetscCall(PetscSectionSetNumFields(coordSection, 1));
2865f34f2dcSJed Brown   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
2875f34f2dcSJed Brown   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
2885f34f2dcSJed Brown   for (v = numCells; v < numCells + numVertices; ++v) {
2895f34f2dcSJed Brown     PetscCall(PetscSectionSetDof(coordSection, v, dim));
2905f34f2dcSJed Brown     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
2915f34f2dcSJed Brown   }
2925f34f2dcSJed Brown   PetscCall(PetscSectionSetUp(coordSection));
2935f34f2dcSJed Brown 
2945f34f2dcSJed Brown   PetscCall(DMCreateLocalVector(cdm, &coordinates));
2955f34f2dcSJed Brown   PetscCall(VecGetArray(coordinates, &coords));
2965f34f2dcSJed Brown   if (rank == 0) {
2975f34f2dcSJed Brown     PetscInt off = 0;
2985f34f2dcSJed Brown     float   *x[3];
2995f34f2dcSJed Brown     int      z, d;
3005f34f2dcSJed Brown 
3015f34f2dcSJed Brown     PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2]));
3025f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
3035f34f2dcSJed Brown       CGNS_ENUMT(DataType_t) datatype;
3045f34f2dcSJed Brown       cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
3055f34f2dcSJed Brown       cgsize_t range_min[3] = {1, 1, 1};
3065f34f2dcSJed Brown       cgsize_t range_max[3] = {1, 1, 1};
3075f34f2dcSJed Brown       int      ngrids, ncoords;
3085f34f2dcSJed Brown 
3095f34f2dcSJed Brown       PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
3105f34f2dcSJed Brown       range_max[0] = sizes[0];
3115f34f2dcSJed Brown       PetscCallCGNS(cg_ngrids(cgid, 1, z, &ngrids));
3125f34f2dcSJed Brown       PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids);
3135f34f2dcSJed Brown       PetscCallCGNS(cg_ncoords(cgid, 1, z, &ncoords));
3145f34f2dcSJed Brown       PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords);
3155f34f2dcSJed Brown       for (d = 0; d < coordDim; ++d) {
3165f34f2dcSJed Brown         PetscCallCGNS(cg_coord_info(cgid, 1, z, 1 + d, &datatype, buffer));
3175f34f2dcSJed Brown         PetscCallCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]));
3185f34f2dcSJed Brown       }
3195f34f2dcSJed Brown       if (coordDim >= 1) {
3205f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v];
3215f34f2dcSJed Brown       }
3225f34f2dcSJed Brown       if (coordDim >= 2) {
3235f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v];
3245f34f2dcSJed Brown       }
3255f34f2dcSJed Brown       if (coordDim >= 3) {
3265f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v];
3275f34f2dcSJed Brown       }
3285f34f2dcSJed Brown       off += sizes[0];
3295f34f2dcSJed Brown     }
3305f34f2dcSJed Brown     PetscCall(PetscFree3(x[0], x[1], x[2]));
3315f34f2dcSJed Brown   }
3325f34f2dcSJed Brown   PetscCall(VecRestoreArray(coordinates, &coords));
3335f34f2dcSJed Brown 
3345f34f2dcSJed Brown   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
3355f34f2dcSJed Brown   PetscCall(VecSetBlockSize(coordinates, coordDim));
3365f34f2dcSJed Brown   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
3375f34f2dcSJed Brown   PetscCall(VecDestroy(&coordinates));
3385f34f2dcSJed Brown 
3395f34f2dcSJed Brown   /* Read boundary conditions */
3405f34f2dcSJed Brown   PetscCall(DMGetNumLabels(*dm, &labelIdRange[0]));
3415f34f2dcSJed Brown   if (rank == 0) {
3425f34f2dcSJed Brown     CGNS_ENUMT(BCType_t) bctype;
3435f34f2dcSJed Brown     CGNS_ENUMT(DataType_t) datatype;
3445f34f2dcSJed Brown     CGNS_ENUMT(PointSetType_t) pointtype;
3455f34f2dcSJed Brown     cgsize_t  *points;
3465f34f2dcSJed Brown     PetscReal *normals;
3475f34f2dcSJed Brown     int        normal[3];
3485f34f2dcSJed Brown     char      *bcname = buffer;
3495f34f2dcSJed Brown     cgsize_t   npoints, nnormals;
3505f34f2dcSJed Brown     int        z, nbc, bc, c, ndatasets;
3515f34f2dcSJed Brown 
3525f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
3535f34f2dcSJed Brown       PetscCallCGNS(cg_nbocos(cgid, 1, z, &nbc));
3545f34f2dcSJed Brown       for (bc = 1; bc <= nbc; ++bc) {
3555f34f2dcSJed Brown         PetscCallCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets));
3565f34f2dcSJed Brown         PetscCall(DMCreateLabel(*dm, bcname));
3575f34f2dcSJed Brown         PetscCall(DMGetLabel(*dm, bcname, &label));
3585f34f2dcSJed Brown         PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals));
3595f34f2dcSJed Brown         PetscCallCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *)normals));
3605f34f2dcSJed Brown         if (pointtype == CGNS_ENUMV(ElementRange)) {
3615f34f2dcSJed Brown           /* Range of cells: assuming half-open interval since the documentation sucks */
36248a46eb9SPierre Jolivet           for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
3635f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(ElementList)) {
3645f34f2dcSJed Brown           /* List of cells */
36548a46eb9SPierre Jolivet           for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
3665f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(PointRange)) {
3675f34f2dcSJed Brown           CGNS_ENUMT(GridLocation_t) gridloc;
3685f34f2dcSJed Brown 
3695f34f2dcSJed Brown           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
3705f34f2dcSJed Brown           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
3715f34f2dcSJed Brown           PetscCallCGNS(cg_gridlocation_read(&gridloc));
3725f34f2dcSJed Brown           /* Range of points: assuming half-open interval since the documentation sucks */
3735f34f2dcSJed Brown           for (c = points[0]; c < points[1]; ++c) {
3745f34f2dcSJed Brown             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1));
3755f34f2dcSJed Brown             else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1));
3765f34f2dcSJed Brown           }
3775f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(PointList)) {
3785f34f2dcSJed Brown           CGNS_ENUMT(GridLocation_t) gridloc;
3795f34f2dcSJed Brown 
3805f34f2dcSJed Brown           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
3815f34f2dcSJed Brown           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
3825f34f2dcSJed Brown           PetscCallCGNS(cg_gridlocation_read(&gridloc));
3835f34f2dcSJed Brown           for (c = 0; c < npoints; ++c) {
3845f34f2dcSJed Brown             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1));
3855f34f2dcSJed Brown             else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1));
3865f34f2dcSJed Brown           }
3875f34f2dcSJed Brown         } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype);
3885f34f2dcSJed Brown         PetscCall(PetscFree2(points, normals));
3895f34f2dcSJed Brown       }
3905f34f2dcSJed Brown     }
3915f34f2dcSJed Brown     PetscCall(PetscFree2(cellStart, vertStart));
3925f34f2dcSJed Brown   }
3935f34f2dcSJed Brown   PetscCall(DMGetNumLabels(*dm, &labelIdRange[1]));
3945f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm));
3955f34f2dcSJed Brown 
3965f34f2dcSJed Brown   /* Create BC labels at all processes */
3975f34f2dcSJed Brown   for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
3985f34f2dcSJed Brown     char       *labelName = buffer;
3995f34f2dcSJed Brown     size_t      len       = sizeof(buffer);
4005f34f2dcSJed Brown     const char *locName;
4015f34f2dcSJed Brown 
4025f34f2dcSJed Brown     if (rank == 0) {
4035f34f2dcSJed Brown       PetscCall(DMGetLabelByNum(*dm, labelId, &label));
4045f34f2dcSJed Brown       PetscCall(PetscObjectGetName((PetscObject)label, &locName));
4055f34f2dcSJed Brown       PetscCall(PetscStrncpy(labelName, locName, len));
4065f34f2dcSJed Brown     }
4075f34f2dcSJed Brown     PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
4085f34f2dcSJed Brown     PetscCallMPI(DMCreateLabel(*dm, labelName));
4095f34f2dcSJed Brown   }
4103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4115f34f2dcSJed Brown }
4125f34f2dcSJed Brown 
4135f34f2dcSJed Brown // Permute plex closure ordering to CGNS
414d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm)
415d71ae5a4SJacob Faibussowitsch {
4165f34f2dcSJed Brown   // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unst_example
4175f34f2dcSJed Brown   static const int bar_2[2]   = {0, 1};
4185f34f2dcSJed Brown   static const int bar_3[3]   = {1, 2, 0};
4195f34f2dcSJed Brown   static const int bar_4[4]   = {2, 3, 0, 1};
4205f34f2dcSJed Brown   static const int bar_5[5]   = {3, 4, 0, 1, 2};
4215f34f2dcSJed Brown   static const int tri_3[3]   = {0, 1, 2};
4225f34f2dcSJed Brown   static const int tri_6[6]   = {3, 4, 5, 0, 1, 2};
4235f34f2dcSJed Brown   static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0};
4245f34f2dcSJed Brown   static const int quad_4[4]  = {0, 1, 2, 3};
4255f34f2dcSJed Brown   static const int quad_9[9]  = {
4265f34f2dcSJed Brown     5, 6, 7, 8, // vertices
4275f34f2dcSJed Brown     1, 2, 3, 4, // edges
4285f34f2dcSJed Brown     0,          // center
4295f34f2dcSJed Brown   };
4305f34f2dcSJed Brown   static const int quad_16[] = {
4315f34f2dcSJed Brown     12, 13, 14, 15,               // vertices
4325f34f2dcSJed Brown     4,  5,  6,  7,  8, 9, 10, 11, // edges
4335f34f2dcSJed Brown     0,  1,  3,  2,                // centers
4345f34f2dcSJed Brown   };
4355f34f2dcSJed Brown   static const int quad_25[] = {
4365f34f2dcSJed Brown     21, 22, 23, 24,                                 // vertices
4375f34f2dcSJed Brown     9,  10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges
4385f34f2dcSJed Brown     0,  1,  2,  5,  8,  7,  6,  3,  4,              // centers
4395f34f2dcSJed Brown   };
4405f34f2dcSJed Brown   static const int tetra_4[4]   = {0, 2, 1, 3};
4415f34f2dcSJed Brown   static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4};
4425f34f2dcSJed Brown   static const int tetra_20[20] = {
4435f34f2dcSJed Brown     16, 18, 17, 19,         // vertices
4445f34f2dcSJed Brown     9,  8,  7,  6,  5,  4,  // bottom edges
4455f34f2dcSJed Brown     10, 11, 14, 15, 13, 12, // side edges
4465f34f2dcSJed Brown     0,  2,  3,  1,          // faces
4475f34f2dcSJed Brown   };
4485f34f2dcSJed Brown   static const int hexa_8[8]   = {0, 3, 2, 1, 4, 5, 6, 7};
4495f34f2dcSJed Brown   static const int hexa_27[27] = {
4505f34f2dcSJed Brown     19, 22, 21, 20, 23, 24, 25, 26, // vertices
4515f34f2dcSJed Brown     10, 9,  8,  7,                  // bottom edges
4525f34f2dcSJed Brown     16, 15, 18, 17,                 // mid edges
4535f34f2dcSJed Brown     11, 12, 13, 14,                 // top edges
4545f34f2dcSJed Brown     1,  3,  5,  4,  6,  2,          // faces
4555f34f2dcSJed Brown     0,                              // center
4565f34f2dcSJed Brown   };
4579371c9d4SSatish Balay   static const int hexa_64[64] = {
4589371c9d4SSatish Balay     // debug with $PETSC_ARCH/tests/dm/impls/plex/tests/ex49 -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 1,1,1 -dm_coord_petscspace_degree 3
4595f34f2dcSJed Brown     56, 59, 58, 57, 60, 61, 62, 63, // vertices
4605f34f2dcSJed Brown     39, 38, 37, 36, 35, 34, 33, 32, // bottom edges
4615f34f2dcSJed Brown     51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24
4625f34f2dcSJed Brown     40, 41, 42, 43, 44, 45, 46, 47, // top edges
4635f34f2dcSJed Brown     8,  10, 11, 9,                  // z-minus face
4645f34f2dcSJed Brown     16, 17, 19, 18,                 // y-minus face
4655f34f2dcSJed Brown     24, 25, 27, 26,                 // x-plus face
4665f34f2dcSJed Brown     20, 21, 23, 22,                 // y-plus face
4675f34f2dcSJed Brown     30, 28, 29, 31,                 // x-minus face
4685f34f2dcSJed Brown     12, 13, 15, 14,                 // z-plus face
4695f34f2dcSJed Brown     0,  1,  3,  2,  4,  5,  7,  6,  // center
4705f34f2dcSJed Brown   };
4715f34f2dcSJed Brown 
4725f34f2dcSJed Brown   PetscFunctionBegin;
4735f34f2dcSJed Brown   *element_type = CGNS_ENUMV(ElementTypeNull);
4745f34f2dcSJed Brown   *perm         = NULL;
4755f34f2dcSJed Brown   switch (cell_type) {
4765f34f2dcSJed Brown   case DM_POLYTOPE_SEGMENT:
4775f34f2dcSJed Brown     switch (closure_size) {
478d71ae5a4SJacob Faibussowitsch     case 2:
479d71ae5a4SJacob Faibussowitsch       *element_type = CGNS_ENUMV(BAR_2);
480d71ae5a4SJacob Faibussowitsch       *perm         = bar_2;
481d71ae5a4SJacob Faibussowitsch     case 3:
482d71ae5a4SJacob Faibussowitsch       *element_type = CGNS_ENUMV(BAR_3);
483d71ae5a4SJacob Faibussowitsch       *perm         = bar_3;
4845f34f2dcSJed Brown     case 4:
4855f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_4);
4865f34f2dcSJed Brown       *perm         = bar_4;
4875f34f2dcSJed Brown       break;
4885f34f2dcSJed Brown     case 5:
4895f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_5);
4905f34f2dcSJed Brown       *perm         = bar_5;
4915f34f2dcSJed Brown       break;
492d71ae5a4SJacob Faibussowitsch     default:
493d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
4945f34f2dcSJed Brown     }
4955f34f2dcSJed Brown     break;
4965f34f2dcSJed Brown   case DM_POLYTOPE_TRIANGLE:
4975f34f2dcSJed Brown     switch (closure_size) {
4985f34f2dcSJed Brown     case 3:
4995f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_3);
5005f34f2dcSJed Brown       *perm         = tri_3;
5015f34f2dcSJed Brown       break;
5025f34f2dcSJed Brown     case 6:
5035f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_6);
5045f34f2dcSJed Brown       *perm         = tri_6;
5055f34f2dcSJed Brown       break;
5065f34f2dcSJed Brown     case 10:
5075f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_10);
5085f34f2dcSJed Brown       *perm         = tri_10;
5095f34f2dcSJed Brown       break;
510d71ae5a4SJacob Faibussowitsch     default:
511d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5125f34f2dcSJed Brown     }
5135f34f2dcSJed Brown     break;
5145f34f2dcSJed Brown   case DM_POLYTOPE_QUADRILATERAL:
5155f34f2dcSJed Brown     switch (closure_size) {
5165f34f2dcSJed Brown     case 4:
5175f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_4);
5185f34f2dcSJed Brown       *perm         = quad_4;
5195f34f2dcSJed Brown       break;
5205f34f2dcSJed Brown     case 9:
5215f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_9);
5225f34f2dcSJed Brown       *perm         = quad_9;
5235f34f2dcSJed Brown       break;
5245f34f2dcSJed Brown     case 16:
5255f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_16);
5265f34f2dcSJed Brown       *perm         = quad_16;
5275f34f2dcSJed Brown       break;
5285f34f2dcSJed Brown     case 25:
5295f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_25);
5305f34f2dcSJed Brown       *perm         = quad_25;
5315f34f2dcSJed Brown       break;
532d71ae5a4SJacob Faibussowitsch     default:
533d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5345f34f2dcSJed Brown     }
5355f34f2dcSJed Brown     break;
5365f34f2dcSJed Brown   case DM_POLYTOPE_TETRAHEDRON:
5375f34f2dcSJed Brown     switch (closure_size) {
5385f34f2dcSJed Brown     case 4:
5395f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TETRA_4);
5405f34f2dcSJed Brown       *perm         = tetra_4;
5415f34f2dcSJed Brown       break;
5425f34f2dcSJed Brown     case 10:
5435f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TETRA_10);
5445f34f2dcSJed Brown       *perm         = tetra_10;
5455f34f2dcSJed Brown       break;
5465f34f2dcSJed Brown     case 20:
5475f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TETRA_20);
5485f34f2dcSJed Brown       *perm         = tetra_20;
5495f34f2dcSJed Brown       break;
550d71ae5a4SJacob Faibussowitsch     default:
551d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5525f34f2dcSJed Brown     }
5535f34f2dcSJed Brown     break;
5545f34f2dcSJed Brown   case DM_POLYTOPE_HEXAHEDRON:
5555f34f2dcSJed Brown     switch (closure_size) {
5565f34f2dcSJed Brown     case 8:
5575f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_8);
5585f34f2dcSJed Brown       *perm         = hexa_8;
5595f34f2dcSJed Brown       break;
5605f34f2dcSJed Brown     case 27:
5615f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_27);
5625f34f2dcSJed Brown       *perm         = hexa_27;
5635f34f2dcSJed Brown       break;
5645f34f2dcSJed Brown     case 64:
5655f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_64);
5665f34f2dcSJed Brown       *perm         = hexa_64;
5675f34f2dcSJed Brown       break;
568d71ae5a4SJacob Faibussowitsch     default:
569d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5705f34f2dcSJed Brown     }
5715f34f2dcSJed Brown     break;
572d71ae5a4SJacob Faibussowitsch   default:
573d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5745f34f2dcSJed Brown   }
5753ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5765f34f2dcSJed Brown }
5775f34f2dcSJed Brown 
5785f34f2dcSJed Brown // node_l2g must be freed
579d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g)
580d71ae5a4SJacob Faibussowitsch {
5815f34f2dcSJed Brown   PetscSection    local_section;
5825f34f2dcSJed Brown   PetscSF         point_sf;
5835f34f2dcSJed Brown   PetscInt        pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes;
5845f34f2dcSJed Brown   PetscMPIInt     comm_size;
5855f34f2dcSJed Brown   const PetscInt *ilocal, field = 0;
5865f34f2dcSJed Brown 
5875f34f2dcSJed Brown   PetscFunctionBegin;
5885f34f2dcSJed Brown   *num_local_nodes  = -1;
5895f34f2dcSJed Brown   *num_global_nodes = -1;
5905f34f2dcSJed Brown   *nStart           = -1;
5915f34f2dcSJed Brown   *nEnd             = -1;
5925f34f2dcSJed Brown   *node_l2g         = NULL;
5935f34f2dcSJed Brown 
5945f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size));
5955f34f2dcSJed Brown   PetscCall(DMGetLocalSection(dm, &local_section));
5965f34f2dcSJed Brown   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5975f34f2dcSJed Brown   PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd));
5985f34f2dcSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf));
5995f34f2dcSJed Brown   if (comm_size == 1) nleaves = 0;
6005f34f2dcSJed Brown   else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL));
6015f34f2dcSJed Brown   PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp));
6025f34f2dcSJed Brown 
6035f34f2dcSJed Brown   PetscInt local_node = 0, owned_node = 0, owned_start = 0;
6045f34f2dcSJed Brown   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
6055f34f2dcSJed Brown     PetscInt dof;
6065f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
6075f34f2dcSJed Brown     PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp);
6085f34f2dcSJed Brown     local_node += dof / ncomp;
6095f34f2dcSJed Brown     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
6105f34f2dcSJed Brown       leaf++;
6115f34f2dcSJed Brown     } else {
6125f34f2dcSJed Brown       owned_node += dof / ncomp;
6135f34f2dcSJed Brown     }
6145f34f2dcSJed Brown   }
6155f34f2dcSJed Brown   PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
6165f34f2dcSJed Brown   PetscCall(PetscMalloc1(pEnd - pStart, &points));
6175f34f2dcSJed Brown   owned_node = 0;
6185f34f2dcSJed Brown   for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) {
6195f34f2dcSJed Brown     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
6205f34f2dcSJed Brown       points[p - pStart] = -1;
6215f34f2dcSJed Brown       leaf++;
6225f34f2dcSJed Brown       continue;
6235f34f2dcSJed Brown     }
6245f34f2dcSJed Brown     PetscInt dof, offset;
6255f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
6265f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
6275f34f2dcSJed Brown     PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp);
6285f34f2dcSJed Brown     points[p - pStart] = owned_start + owned_node;
6295f34f2dcSJed Brown     owned_node += dof / ncomp;
6305f34f2dcSJed Brown   }
6315f34f2dcSJed Brown   if (comm_size > 1) {
6325f34f2dcSJed Brown     PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE));
6335f34f2dcSJed Brown     PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE));
6345f34f2dcSJed Brown   }
6355f34f2dcSJed Brown 
6365f34f2dcSJed Brown   // Set up global indices for each local node
6375f34f2dcSJed Brown   PetscCall(PetscMalloc1(local_node, &nodes));
6385f34f2dcSJed Brown   for (PetscInt p = spStart; p < spEnd; p++) {
6395f34f2dcSJed Brown     PetscInt dof, offset;
6405f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
6415f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
6425f34f2dcSJed Brown     for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n;
6435f34f2dcSJed Brown   }
6445f34f2dcSJed Brown   PetscCall(PetscFree(points));
6455f34f2dcSJed Brown   *num_local_nodes = local_node;
6465f34f2dcSJed Brown   *nStart          = owned_start;
6475f34f2dcSJed Brown   *nEnd            = owned_start + owned_node;
648712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
6495f34f2dcSJed Brown   *node_l2g = nodes;
6503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6515f34f2dcSJed Brown }
6525f34f2dcSJed Brown 
653d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer)
654d71ae5a4SJacob Faibussowitsch {
6555f34f2dcSJed Brown   PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data;
6565f34f2dcSJed Brown   PetscInt          topo_dim, coord_dim, num_global_elems;
6575f34f2dcSJed Brown   PetscInt          cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd;
6585f34f2dcSJed Brown   const PetscInt   *node_l2g;
6595f34f2dcSJed Brown   Vec               coord;
660b85bf5edSJed Brown   DM                colloc_dm, cdm;
6615f34f2dcSJed Brown   PetscMPIInt       size;
6625f34f2dcSJed Brown   const char       *dm_name;
6635f34f2dcSJed Brown   int               base, zone;
6645f34f2dcSJed Brown   cgsize_t          isize[3];
6655f34f2dcSJed Brown 
6665f34f2dcSJed Brown   PetscFunctionBegin;
66725760affSJed Brown   if (!cgv->file_num) {
66825760affSJed Brown     PetscInt time_step;
66925760affSJed Brown     PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL));
67025760affSJed Brown     PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step));
67125760affSJed Brown   }
6725f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
6735f34f2dcSJed Brown   PetscCall(DMGetDimension(dm, &topo_dim));
6745f34f2dcSJed Brown   PetscCall(DMGetCoordinateDim(dm, &coord_dim));
6755f34f2dcSJed Brown   PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name));
6765f34f2dcSJed Brown   PetscCallCGNS(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base));
6775f34f2dcSJed Brown   PetscCallCGNS(cg_goto(cgv->file_num, base, NULL));
6785f34f2dcSJed Brown   PetscCallCGNS(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)));
6795f34f2dcSJed Brown 
680b85bf5edSJed Brown   {
681b85bf5edSJed Brown     PetscFE        fe, fe_coord;
6829bb8f83bSJed Brown     PetscClassId   ds_id;
683b85bf5edSJed Brown     PetscDualSpace dual_space, dual_space_coord;
6842d1fc720SJed Brown     PetscInt       num_fields, field_order = -1, field_order_coord;
685b85bf5edSJed Brown     PetscBool      is_simplex;
6865b8b2c14SJed Brown     PetscCall(DMGetNumFields(dm, &num_fields));
6879bb8f83bSJed Brown     if (num_fields > 0) {
6889bb8f83bSJed Brown       PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe));
6899bb8f83bSJed Brown       PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id));
6902d1fc720SJed Brown       if (ds_id != PETSCFE_CLASSID) {
6912d1fc720SJed Brown         fe = NULL;
6922d1fc720SJed Brown         if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter
6932d1fc720SJed Brown         else field_order = 1;                           // assume vertex-based linear elements
6942d1fc720SJed Brown       }
6959bb8f83bSJed Brown     } else fe = NULL;
696b85bf5edSJed Brown     if (fe) {
697b85bf5edSJed Brown       PetscCall(PetscFEGetDualSpace(fe, &dual_space));
698b85bf5edSJed Brown       PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order));
6992d1fc720SJed Brown     }
7005f34f2dcSJed Brown     PetscCall(DMGetCoordinateDM(dm, &cdm));
701b85bf5edSJed Brown     PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord));
7029bb8f83bSJed Brown     {
7039bb8f83bSJed Brown       PetscClassId id;
7049bb8f83bSJed Brown       PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id));
7059bb8f83bSJed Brown       if (id != PETSCFE_CLASSID) fe_coord = NULL;
7069bb8f83bSJed Brown     }
707b85bf5edSJed Brown     if (fe_coord) {
708b85bf5edSJed Brown       PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord));
709b85bf5edSJed Brown       PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord));
710b85bf5edSJed Brown     } else field_order_coord = 1;
7112d1fc720SJed Brown     if (field_order > 0 && field_order != field_order_coord) {
712b85bf5edSJed Brown       PetscInt quadrature_order = field_order;
713b85bf5edSJed Brown       PetscCall(DMClone(dm, &colloc_dm));
7147d4eb7abSJed Brown       { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied
715*1fca310dSJames Wright         const PetscSF *face_sfs;
716*1fca310dSJames Wright         PetscInt       num_face_sfs;
717*1fca310dSJames Wright         PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs));
718*1fca310dSJames Wright         PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs));
719*1fca310dSJames Wright         if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal;
7207d4eb7abSJed Brown       }
721b85bf5edSJed Brown       PetscCall(DMPlexIsSimplex(dm, &is_simplex));
722b85bf5edSJed Brown       PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe));
723e44f6aebSMatthew G. Knepley       PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_TRUE));
724b85bf5edSJed Brown       PetscCall(PetscFEDestroy(&fe));
725b85bf5edSJed Brown     } else {
726b85bf5edSJed Brown       PetscCall(PetscObjectReference((PetscObject)dm));
727b85bf5edSJed Brown       colloc_dm = dm;
728b85bf5edSJed Brown     }
729b85bf5edSJed Brown   }
730b85bf5edSJed Brown   PetscCall(DMGetCoordinateDM(colloc_dm, &cdm));
7315f34f2dcSJed Brown   PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g));
732b85bf5edSJed Brown   PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord));
7335f34f2dcSJed Brown   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
7345f34f2dcSJed Brown   num_global_elems = cEnd - cStart;
735712fec58SPierre Jolivet   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
7365f34f2dcSJed Brown   isize[0] = num_global_nodes;
7375f34f2dcSJed Brown   isize[1] = num_global_elems;
7385f34f2dcSJed Brown   isize[2] = 0;
7395f34f2dcSJed Brown   PetscCallCGNS(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone));
7405f34f2dcSJed Brown 
7415f34f2dcSJed Brown   {
7425f34f2dcSJed Brown     const PetscScalar *X;
7435f34f2dcSJed Brown     PetscScalar       *x;
7445f34f2dcSJed Brown     int                coord_ids[3];
7455f34f2dcSJed Brown 
7465f34f2dcSJed Brown     PetscCall(VecGetArrayRead(coord, &X));
7475f34f2dcSJed Brown     for (PetscInt d = 0; d < coord_dim; d++) {
7485f34f2dcSJed Brown       const double exponents[] = {0, 1, 0, 0, 0};
7495f34f2dcSJed Brown       char         coord_name[64];
7505f34f2dcSJed Brown       PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d));
7515f34f2dcSJed Brown       PetscCallCGNS(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]));
7525f34f2dcSJed Brown       PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL));
7535f34f2dcSJed Brown       PetscCallCGNS(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents));
7545f34f2dcSJed Brown     }
7555f34f2dcSJed Brown 
7565f34f2dcSJed Brown     DMPolytopeType cell_type;
7575f34f2dcSJed Brown     int            section;
7585f34f2dcSJed Brown     cgsize_t       e_owned, e_global, e_start, *conn = NULL;
7595f34f2dcSJed Brown     const int     *perm;
760e853fb4cSJames Wright     CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull);
7615f34f2dcSJed Brown     {
7625f34f2dcSJed Brown       PetscCall(PetscMalloc1(nEnd - nStart, &x));
7635f34f2dcSJed Brown       for (PetscInt d = 0; d < coord_dim; d++) {
7645f34f2dcSJed Brown         for (PetscInt n = 0; n < num_local_nodes; n++) {
7655f34f2dcSJed Brown           PetscInt gn = node_l2g[n];
7665f34f2dcSJed Brown           if (gn < nStart || nEnd <= gn) continue;
7675f34f2dcSJed Brown           x[gn - nStart] = X[n * coord_dim + d];
7685f34f2dcSJed Brown         }
7695f34f2dcSJed Brown         // CGNS nodes use 1-based indexing
7705f34f2dcSJed Brown         cgsize_t start = nStart + 1, end = nEnd;
7715f34f2dcSJed Brown         PetscCallCGNS(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x));
7725f34f2dcSJed Brown       }
7735f34f2dcSJed Brown       PetscCall(PetscFree(x));
7745f34f2dcSJed Brown       PetscCall(VecRestoreArrayRead(coord, &X));
7755f34f2dcSJed Brown     }
7765f34f2dcSJed Brown 
7775f34f2dcSJed Brown     PetscCall(DMPlexGetCellType(dm, cStart, &cell_type));
7785f34f2dcSJed Brown     for (PetscInt i = cStart, c = 0; i < cEnd; i++) {
7795f34f2dcSJed Brown       PetscInt closure_dof, *closure_indices, elem_size;
7805f34f2dcSJed Brown       PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
7815f34f2dcSJed Brown       elem_size = closure_dof / coord_dim;
7825f34f2dcSJed Brown       if (!conn) PetscCall(PetscMalloc1((cEnd - cStart) * elem_size, &conn));
7835f34f2dcSJed Brown       PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm));
784ad540459SPierre Jolivet       for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1;
7855f34f2dcSJed Brown       PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
7865f34f2dcSJed Brown     }
7875f34f2dcSJed Brown     e_owned = cEnd - cStart;
788712fec58SPierre Jolivet     PetscCall(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm)));
7895f34f2dcSJed Brown     PetscCheck(e_global == num_global_elems, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected number of elements %" PetscInt64_FMT "vs %" PetscInt_FMT, e_global, num_global_elems);
7905f34f2dcSJed Brown     e_start = 0;
7915f34f2dcSJed Brown     PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm)));
7925f34f2dcSJed Brown     PetscCallCGNS(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, &section));
7935f34f2dcSJed Brown     PetscCallCGNS(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn));
7945f34f2dcSJed Brown     PetscCall(PetscFree(conn));
7955f34f2dcSJed Brown 
7965f34f2dcSJed Brown     cgv->base            = base;
7975f34f2dcSJed Brown     cgv->zone            = zone;
7985f34f2dcSJed Brown     cgv->node_l2g        = node_l2g;
7995f34f2dcSJed Brown     cgv->num_local_nodes = num_local_nodes;
8005f34f2dcSJed Brown     cgv->nStart          = nStart;
8015f34f2dcSJed Brown     cgv->nEnd            = nEnd;
8029bb8f83bSJed Brown     cgv->eStart          = e_start;
8039bb8f83bSJed Brown     cgv->eEnd            = e_start + e_owned;
8045f34f2dcSJed Brown     if (1) {
8055f34f2dcSJed Brown       PetscMPIInt rank;
8065f34f2dcSJed Brown       int        *efield;
8075f34f2dcSJed Brown       int         sol, field;
8085f34f2dcSJed Brown       DMLabel     label;
8095f34f2dcSJed Brown       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
8105f34f2dcSJed Brown       PetscCall(PetscMalloc1(e_owned, &efield));
8115f34f2dcSJed Brown       for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank;
8125f34f2dcSJed Brown       PetscCallCGNS(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol));
8135f34f2dcSJed Brown       PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field));
8145f34f2dcSJed Brown       cgsize_t start = e_start + 1, end = e_start + e_owned;
8155f34f2dcSJed Brown       PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
8165f34f2dcSJed Brown       PetscCall(DMGetLabel(dm, "Cell Sets", &label));
8175f34f2dcSJed Brown       if (label) {
8185f34f2dcSJed Brown         for (PetscInt c = cStart; c < cEnd; c++) {
8195f34f2dcSJed Brown           PetscInt value;
8205f34f2dcSJed Brown           PetscCall(DMLabelGetValue(label, c, &value));
8215f34f2dcSJed Brown           efield[c - cStart] = value;
8225f34f2dcSJed Brown         }
8235f34f2dcSJed Brown         PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field));
8245f34f2dcSJed Brown         PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
8255f34f2dcSJed Brown       }
8265f34f2dcSJed Brown       PetscCall(PetscFree(efield));
8275f34f2dcSJed Brown     }
8285f34f2dcSJed Brown   }
829b85bf5edSJed Brown   PetscCall(DMDestroy(&colloc_dm));
8303ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8315f34f2dcSJed Brown }
8325f34f2dcSJed Brown 
833d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd)
834d71ae5a4SJacob Faibussowitsch {
8355f34f2dcSJed Brown   PetscFunctionBegin;
8365f34f2dcSJed Brown   switch (pd) {
837d71ae5a4SJacob Faibussowitsch   case PETSC_FLOAT:
838d71ae5a4SJacob Faibussowitsch     *cd = CGNS_ENUMV(RealSingle);
839d71ae5a4SJacob Faibussowitsch     break;
840d71ae5a4SJacob Faibussowitsch   case PETSC_DOUBLE:
841d71ae5a4SJacob Faibussowitsch     *cd = CGNS_ENUMV(RealDouble);
842d71ae5a4SJacob Faibussowitsch     break;
843d71ae5a4SJacob Faibussowitsch   case PETSC_COMPLEX:
844d71ae5a4SJacob Faibussowitsch     *cd = CGNS_ENUMV(ComplexDouble);
845d71ae5a4SJacob Faibussowitsch     break;
846d71ae5a4SJacob Faibussowitsch   default:
847d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]);
8485f34f2dcSJed Brown   }
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8505f34f2dcSJed Brown }
8515f34f2dcSJed Brown 
852d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer)
853d71ae5a4SJacob Faibussowitsch {
8545f34f2dcSJed Brown   PetscViewer_CGNS  *cgv = (PetscViewer_CGNS *)viewer->data;
8555f34f2dcSJed Brown   DM                 dm;
8565f34f2dcSJed Brown   PetscSection       section;
8579bb8f83bSJed Brown   PetscInt           time_step, num_fields, pStart, pEnd, cStart, cEnd;
8585f34f2dcSJed Brown   PetscReal          time, *time_slot;
8599812b6beSJed Brown   size_t            *step_slot;
8605f34f2dcSJed Brown   const PetscScalar *v;
8615f34f2dcSJed Brown   char               solution_name[PETSC_MAX_PATH_LEN];
8625f34f2dcSJed Brown   int                sol;
8635f34f2dcSJed Brown 
8645f34f2dcSJed Brown   PetscFunctionBegin;
8655f34f2dcSJed Brown   PetscCall(VecGetDM(V, &dm));
8665f34f2dcSJed Brown   if (!cgv->node_l2g) PetscCall(DMView(dm, viewer));
8679bb8f83bSJed Brown   if (!cgv->nodal_field) PetscCall(PetscMalloc1(PetscMax(cgv->nEnd - cgv->nStart, cgv->eEnd - cgv->eStart), &cgv->nodal_field));
8685f34f2dcSJed Brown   if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times));
8699812b6beSJed Brown   if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps));
8705f34f2dcSJed Brown 
8715f34f2dcSJed Brown   PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time));
8729371c9d4SSatish Balay   if (time_step < 0) {
8739371c9d4SSatish Balay     time_step = 0;
8749371c9d4SSatish Balay     time      = 0.;
8759371c9d4SSatish Balay   }
8765f34f2dcSJed Brown   PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot));
8775f34f2dcSJed Brown   *time_slot = time;
8789812b6beSJed Brown   PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot));
8799812b6beSJed Brown   *step_slot = time_step;
8805f34f2dcSJed Brown   PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step));
8815f34f2dcSJed Brown   PetscCall(DMGetLocalSection(dm, &section));
8829bb8f83bSJed Brown   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
88300a86070SJed Brown   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
8849bb8f83bSJed Brown   CGNS_ENUMT(GridLocation_t) grid_loc = CGNS_ENUMV(Vertex);
8859bb8f83bSJed Brown   if (cStart == pStart && cEnd == pEnd) grid_loc = CGNS_ENUMV(CellCenter);
8869bb8f83bSJed Brown   PetscCallCGNS(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, grid_loc, &sol));
8875f34f2dcSJed Brown   PetscCall(VecGetArrayRead(V, &v));
88800a86070SJed Brown   PetscCall(PetscSectionGetNumFields(section, &num_fields));
88900a86070SJed Brown   for (PetscInt field = 0; field < num_fields; field++) {
89000a86070SJed Brown     PetscInt    ncomp;
8919bb8f83bSJed Brown     const char *field_name;
8929bb8f83bSJed Brown     PetscCall(PetscSectionGetFieldName(section, field, &field_name));
89300a86070SJed Brown     PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp));
8945f34f2dcSJed Brown     for (PetscInt comp = 0; comp < ncomp; comp++) {
8955f34f2dcSJed Brown       int         cgfield;
8965f34f2dcSJed Brown       const char *comp_name;
8979bb8f83bSJed Brown       char        cgns_field_name[32]; // CGNS max field name is 32
8985f34f2dcSJed Brown       CGNS_ENUMT(DataType_t) datatype;
8995f34f2dcSJed Brown       PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name));
9009bb8f83bSJed Brown       if (ncomp == 1 && comp_name[0] == '0' && comp_name[1] == '\0') PetscCall(PetscStrncpy(cgns_field_name, field_name, sizeof cgns_field_name));
9019bb8f83bSJed Brown       else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name));
9029bb8f83bSJed Brown       else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name));
9035f34f2dcSJed Brown       PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype));
9049bb8f83bSJed Brown       PetscCallCGNS(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield));
90500a86070SJed Brown       for (PetscInt p = pStart, n = 0; p < pEnd; p++) {
90600a86070SJed Brown         PetscInt off, dof;
90700a86070SJed Brown         PetscCall(PetscSectionGetFieldDof(section, p, field, &dof));
90800a86070SJed Brown         if (dof == 0) continue;
90900a86070SJed Brown         PetscCall(PetscSectionGetFieldOffset(section, p, field, &off));
91000a86070SJed Brown         for (PetscInt c = comp; c < dof; c += ncomp, n++) {
9119bb8f83bSJed Brown           switch (grid_loc) {
9129bb8f83bSJed Brown           case CGNS_ENUMV(Vertex): {
9135f34f2dcSJed Brown             PetscInt gn = cgv->node_l2g[n];
9145f34f2dcSJed Brown             if (gn < cgv->nStart || cgv->nEnd <= gn) continue;
91500a86070SJed Brown             cgv->nodal_field[gn - cgv->nStart] = v[off + c];
9169bb8f83bSJed Brown           } break;
9179bb8f83bSJed Brown           case CGNS_ENUMV(CellCenter): {
9189bb8f83bSJed Brown             cgv->nodal_field[n] = v[off + c];
9199bb8f83bSJed Brown           } break;
9209bb8f83bSJed Brown           default:
9219bb8f83bSJed Brown             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations");
9229bb8f83bSJed Brown           }
92300a86070SJed Brown         }
9245f34f2dcSJed Brown       }
9255f34f2dcSJed Brown       // CGNS nodes use 1-based indexing
9265f34f2dcSJed Brown       cgsize_t start = cgv->nStart + 1, end = cgv->nEnd;
9279bb8f83bSJed Brown       if (grid_loc == CGNS_ENUMV(CellCenter)) {
9289bb8f83bSJed Brown         start = cgv->eStart + 1;
9299bb8f83bSJed Brown         end   = cgv->eEnd;
9309bb8f83bSJed Brown       }
9315f34f2dcSJed Brown       PetscCallCGNS(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field));
9325f34f2dcSJed Brown     }
93300a86070SJed Brown   }
9345f34f2dcSJed Brown   PetscCall(VecRestoreArrayRead(V, &v));
93525760affSJed Brown   PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer));
9363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9375f34f2dcSJed Brown }
938