xref: /petsc/src/dm/impls/plex/cgns/plexcgns2.c (revision e853fb4c8b5febe3904ee8ab6dd9d8f36768265c)
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 
155f34f2dcSJed Brown PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
165f34f2dcSJed Brown {
175f34f2dcSJed Brown   PetscMPIInt    rank;
185f34f2dcSJed Brown   int cgid = -1;
195f34f2dcSJed Brown 
205f34f2dcSJed Brown   PetscFunctionBegin;
215f34f2dcSJed Brown   PetscValidCharPointer(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));
295f34f2dcSJed Brown   PetscFunctionReturn(0);
305f34f2dcSJed Brown }
315f34f2dcSJed Brown 
325f34f2dcSJed Brown PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm)
335f34f2dcSJed Brown {
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     }
725f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
735f34f2dcSJed Brown       vertStart[z] += numCells;
745f34f2dcSJed Brown     }
755f34f2dcSJed Brown     coordDim = dim;
765f34f2dcSJed Brown   }
775f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm));
785f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
795f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm));
805f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm));
815f34f2dcSJed Brown 
825f34f2dcSJed Brown   PetscCall(PetscObjectSetName((PetscObject) *dm, basename));
835f34f2dcSJed Brown   PetscCall(DMSetDimension(*dm, dim));
845f34f2dcSJed Brown   PetscCall(DMCreateLabel(*dm, "celltype"));
855f34f2dcSJed Brown   PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices));
865f34f2dcSJed Brown 
875f34f2dcSJed Brown   /* Read zone information */
885f34f2dcSJed Brown   if (rank == 0) {
895f34f2dcSJed Brown     int z, c, c_loc;
905f34f2dcSJed Brown 
915f34f2dcSJed Brown     /* Read the cell set connectivity table and build mesh topology
925f34f2dcSJed Brown        CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */
935f34f2dcSJed Brown     /* First set sizes */
945f34f2dcSJed Brown     for (z = 1, c = 0; z <= nzones; ++z) {
955f34f2dcSJed Brown       CGNS_ENUMT(ZoneType_t)    zonetype;
965f34f2dcSJed Brown       int                       nsections;
975f34f2dcSJed Brown       CGNS_ENUMT(ElementType_t) cellType;
985f34f2dcSJed Brown       cgsize_t                  start, end;
995f34f2dcSJed Brown       int                       nbndry, parentFlag;
1005f34f2dcSJed Brown       PetscInt                  numCorners;
1015f34f2dcSJed Brown       DMPolytopeType            ctype;
1025f34f2dcSJed Brown 
1035f34f2dcSJed Brown       PetscCallCGNS(cg_zone_type(cgid, 1, z, &zonetype));
1045f34f2dcSJed Brown       PetscCheck(zonetype != CGNS_ENUMV(Structured),PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS");
1055f34f2dcSJed Brown       PetscCallCGNS(cg_nsections(cgid, 1, z, &nsections));
1065f34f2dcSJed Brown       PetscCheck(nsections <= 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d",nsections);
1075f34f2dcSJed Brown       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
1085f34f2dcSJed Brown       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
1095f34f2dcSJed Brown       if (cellType == CGNS_ENUMV(MIXED)) {
1105f34f2dcSJed Brown         cgsize_t elementDataSize, *elements;
1115f34f2dcSJed Brown         PetscInt off;
1125f34f2dcSJed Brown 
1135f34f2dcSJed Brown         PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
1145f34f2dcSJed Brown         PetscCall(PetscMalloc1(elementDataSize, &elements));
1155f34f2dcSJed Brown         PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
1165f34f2dcSJed Brown         for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) {
1175f34f2dcSJed Brown           switch (elements[off]) {
1185f34f2dcSJed Brown           case CGNS_ENUMV(BAR_2):   numCorners = 2; ctype = DM_POLYTOPE_SEGMENT;       break;
1195f34f2dcSJed Brown           case CGNS_ENUMV(TRI_3):   numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE;      break;
1205f34f2dcSJed Brown           case CGNS_ENUMV(QUAD_4):  numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break;
1215f34f2dcSJed Brown           case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON;   break;
1225f34f2dcSJed Brown           case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM;     break;
1235f34f2dcSJed Brown           case CGNS_ENUMV(HEXA_8):  numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON;    break;
1245f34f2dcSJed Brown           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]);
1255f34f2dcSJed Brown           }
1265f34f2dcSJed Brown           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
1275f34f2dcSJed Brown           PetscCall(DMPlexSetCellType(*dm, c, ctype));
1285f34f2dcSJed Brown           off += numCorners+1;
1295f34f2dcSJed Brown         }
1305f34f2dcSJed Brown         PetscCall(PetscFree(elements));
1315f34f2dcSJed Brown       } else {
1325f34f2dcSJed Brown         switch (cellType) {
1335f34f2dcSJed Brown         case CGNS_ENUMV(BAR_2):   numCorners = 2; ctype = DM_POLYTOPE_SEGMENT;       break;
1345f34f2dcSJed Brown         case CGNS_ENUMV(TRI_3):   numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE;      break;
1355f34f2dcSJed Brown         case CGNS_ENUMV(QUAD_4):  numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break;
1365f34f2dcSJed Brown         case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON;   break;
1375f34f2dcSJed Brown         case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM;     break;
1385f34f2dcSJed Brown         case CGNS_ENUMV(HEXA_8):  numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON;    break;
1395f34f2dcSJed Brown         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
1405f34f2dcSJed Brown         }
1415f34f2dcSJed Brown         for (c_loc = start; c_loc <= end; ++c_loc, ++c) {
1425f34f2dcSJed Brown           PetscCall(DMPlexSetConeSize(*dm, c, numCorners));
1435f34f2dcSJed Brown           PetscCall(DMPlexSetCellType(*dm, c, ctype));
1445f34f2dcSJed Brown         }
1455f34f2dcSJed Brown       }
1465f34f2dcSJed Brown     }
1475f34f2dcSJed Brown     for (v = numCells; v < numCells+numVertices; ++v) {
1485f34f2dcSJed Brown       PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
1495f34f2dcSJed Brown     }
1505f34f2dcSJed Brown   }
1515f34f2dcSJed Brown 
1525f34f2dcSJed Brown   PetscCall(DMSetUp(*dm));
1535f34f2dcSJed Brown 
1545f34f2dcSJed Brown   PetscCall(DMCreateLabel(*dm, "zone"));
1555f34f2dcSJed Brown   if (rank == 0) {
1565f34f2dcSJed Brown     int z, c, c_loc, v_loc;
1575f34f2dcSJed Brown 
1585f34f2dcSJed Brown     PetscCall(DMGetLabel(*dm, "zone", &label));
1595f34f2dcSJed Brown     for (z = 1, c = 0; z <= nzones; ++z) {
1605f34f2dcSJed Brown       CGNS_ENUMT(ElementType_t)   cellType;
1615f34f2dcSJed Brown       cgsize_t                    elementDataSize, *elements, start, end;
1625f34f2dcSJed Brown       int                          nbndry, parentFlag;
1635f34f2dcSJed Brown       PetscInt                    *cone, numc, numCorners, maxCorners = 27;
1645f34f2dcSJed Brown 
1655f34f2dcSJed Brown       PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag));
1665f34f2dcSJed Brown       numc = end - start;
1675f34f2dcSJed Brown       /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */
1685f34f2dcSJed Brown       PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize));
1695f34f2dcSJed Brown       PetscCall(PetscMalloc2(elementDataSize,&elements,maxCorners,&cone));
1705f34f2dcSJed Brown       PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL));
1715f34f2dcSJed Brown       if (cellType == CGNS_ENUMV(MIXED)) {
1725f34f2dcSJed Brown         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1735f34f2dcSJed Brown         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
1745f34f2dcSJed Brown           switch (elements[v]) {
1755f34f2dcSJed Brown           case CGNS_ENUMV(BAR_2):   numCorners = 2; break;
1765f34f2dcSJed Brown           case CGNS_ENUMV(TRI_3):   numCorners = 3; break;
1775f34f2dcSJed Brown           case CGNS_ENUMV(QUAD_4):  numCorners = 4; break;
1785f34f2dcSJed Brown           case CGNS_ENUMV(TETRA_4): numCorners = 4; break;
1795f34f2dcSJed Brown           case CGNS_ENUMV(PENTA_6): numCorners = 6; break;
1805f34f2dcSJed Brown           case CGNS_ENUMV(HEXA_8):  numCorners = 8; break;
1815f34f2dcSJed Brown           default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]);
1825f34f2dcSJed Brown           }
1835f34f2dcSJed Brown           ++v;
1845f34f2dcSJed Brown           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
1855f34f2dcSJed Brown             cone[v_loc] = elements[v]+numCells-1;
1865f34f2dcSJed Brown           }
1875f34f2dcSJed Brown           PetscCall(DMPlexReorderCell(*dm, c, cone));
1885f34f2dcSJed Brown           PetscCall(DMPlexSetCone(*dm, c, cone));
1895f34f2dcSJed Brown           PetscCall(DMLabelSetValue(label, c, z));
1905f34f2dcSJed Brown         }
1915f34f2dcSJed Brown       } else {
1925f34f2dcSJed Brown         switch (cellType) {
1935f34f2dcSJed Brown         case CGNS_ENUMV(BAR_2):   numCorners = 2; break;
1945f34f2dcSJed Brown         case CGNS_ENUMV(TRI_3):   numCorners = 3; break;
1955f34f2dcSJed Brown         case CGNS_ENUMV(QUAD_4):  numCorners = 4; break;
1965f34f2dcSJed Brown         case CGNS_ENUMV(TETRA_4): numCorners = 4; break;
1975f34f2dcSJed Brown         case CGNS_ENUMV(PENTA_6): numCorners = 6; break;
1985f34f2dcSJed Brown         case CGNS_ENUMV(HEXA_8):  numCorners = 8; break;
1995f34f2dcSJed Brown         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType);
2005f34f2dcSJed Brown         }
2015f34f2dcSJed Brown         /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
2025f34f2dcSJed Brown         for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) {
2035f34f2dcSJed Brown           for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) {
2045f34f2dcSJed Brown             cone[v_loc] = elements[v]+numCells-1;
2055f34f2dcSJed Brown           }
2065f34f2dcSJed Brown           PetscCall(DMPlexReorderCell(*dm, c, cone));
2075f34f2dcSJed Brown           PetscCall(DMPlexSetCone(*dm, c, cone));
2085f34f2dcSJed Brown           PetscCall(DMLabelSetValue(label, c, z));
2095f34f2dcSJed Brown         }
2105f34f2dcSJed Brown       }
2115f34f2dcSJed Brown       PetscCall(PetscFree2(elements,cone));
2125f34f2dcSJed Brown     }
2135f34f2dcSJed Brown   }
2145f34f2dcSJed Brown 
2155f34f2dcSJed Brown   PetscCall(DMPlexSymmetrize(*dm));
2165f34f2dcSJed Brown   PetscCall(DMPlexStratify(*dm));
2175f34f2dcSJed Brown   if (interpolate) {
2185f34f2dcSJed Brown     DM idm;
2195f34f2dcSJed Brown 
2205f34f2dcSJed Brown     PetscCall(DMPlexInterpolate(*dm, &idm));
2215f34f2dcSJed Brown     PetscCall(DMDestroy(dm));
2225f34f2dcSJed Brown     *dm  = idm;
2235f34f2dcSJed Brown   }
2245f34f2dcSJed Brown 
2255f34f2dcSJed Brown   /* Read coordinates */
2265f34f2dcSJed Brown   PetscCall(DMSetCoordinateDim(*dm, coordDim));
2275f34f2dcSJed Brown   PetscCall(DMGetCoordinateDM(*dm, &cdm));
2285f34f2dcSJed Brown   PetscCall(DMGetLocalSection(cdm, &coordSection));
2295f34f2dcSJed Brown   PetscCall(PetscSectionSetNumFields(coordSection, 1));
2305f34f2dcSJed Brown   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim));
2315f34f2dcSJed Brown   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
2325f34f2dcSJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
2335f34f2dcSJed Brown     PetscCall(PetscSectionSetDof(coordSection, v, dim));
2345f34f2dcSJed Brown     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim));
2355f34f2dcSJed Brown   }
2365f34f2dcSJed Brown   PetscCall(PetscSectionSetUp(coordSection));
2375f34f2dcSJed Brown 
2385f34f2dcSJed Brown   PetscCall(DMCreateLocalVector(cdm, &coordinates));
2395f34f2dcSJed Brown   PetscCall(VecGetArray(coordinates, &coords));
2405f34f2dcSJed Brown   if (rank == 0) {
2415f34f2dcSJed Brown     PetscInt off = 0;
2425f34f2dcSJed Brown     float   *x[3];
2435f34f2dcSJed Brown     int      z, d;
2445f34f2dcSJed Brown 
2455f34f2dcSJed Brown     PetscCall(PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]));
2465f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
2475f34f2dcSJed Brown       CGNS_ENUMT(DataType_t) datatype;
2485f34f2dcSJed Brown       cgsize_t               sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */
2495f34f2dcSJed Brown       cgsize_t               range_min[3] = {1, 1, 1};
2505f34f2dcSJed Brown       cgsize_t               range_max[3] = {1, 1, 1};
2515f34f2dcSJed Brown       int                    ngrids, ncoords;
2525f34f2dcSJed Brown 
2535f34f2dcSJed Brown       PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes));
2545f34f2dcSJed Brown       range_max[0] = sizes[0];
2555f34f2dcSJed Brown       PetscCallCGNS(cg_ngrids(cgid, 1, z, &ngrids));
2565f34f2dcSJed Brown       PetscCheck(ngrids <= 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d",ngrids);
2575f34f2dcSJed Brown       PetscCallCGNS(cg_ncoords(cgid, 1, z, &ncoords));
2585f34f2dcSJed Brown       PetscCheck(ncoords == coordDim,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d",ncoords);
2595f34f2dcSJed Brown       for (d = 0; d < coordDim; ++d) {
2605f34f2dcSJed Brown         PetscCallCGNS(cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer));
2615f34f2dcSJed Brown         PetscCallCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]));
2625f34f2dcSJed Brown       }
2635f34f2dcSJed Brown       if (coordDim >= 1) {
2645f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v];
2655f34f2dcSJed Brown       }
2665f34f2dcSJed Brown       if (coordDim >= 2) {
2675f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v];
2685f34f2dcSJed Brown       }
2695f34f2dcSJed Brown       if (coordDim >= 3) {
2705f34f2dcSJed Brown         for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v];
2715f34f2dcSJed Brown       }
2725f34f2dcSJed Brown       off += sizes[0];
2735f34f2dcSJed Brown     }
2745f34f2dcSJed Brown     PetscCall(PetscFree3(x[0],x[1],x[2]));
2755f34f2dcSJed Brown   }
2765f34f2dcSJed Brown   PetscCall(VecRestoreArray(coordinates, &coords));
2775f34f2dcSJed Brown 
2785f34f2dcSJed Brown   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
2795f34f2dcSJed Brown   PetscCall(VecSetBlockSize(coordinates, coordDim));
2805f34f2dcSJed Brown   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
2815f34f2dcSJed Brown   PetscCall(VecDestroy(&coordinates));
2825f34f2dcSJed Brown 
2835f34f2dcSJed Brown   /* Read boundary conditions */
2845f34f2dcSJed Brown   PetscCall(DMGetNumLabels(*dm, &labelIdRange[0]));
2855f34f2dcSJed Brown   if (rank == 0) {
2865f34f2dcSJed Brown     CGNS_ENUMT(BCType_t)        bctype;
2875f34f2dcSJed Brown     CGNS_ENUMT(DataType_t)      datatype;
2885f34f2dcSJed Brown     CGNS_ENUMT(PointSetType_t)  pointtype;
2895f34f2dcSJed Brown     cgsize_t                    *points;
2905f34f2dcSJed Brown     PetscReal                   *normals;
2915f34f2dcSJed Brown     int                         normal[3];
2925f34f2dcSJed Brown     char                        *bcname = buffer;
2935f34f2dcSJed Brown     cgsize_t                    npoints, nnormals;
2945f34f2dcSJed Brown     int                         z, nbc, bc, c, ndatasets;
2955f34f2dcSJed Brown 
2965f34f2dcSJed Brown     for (z = 1; z <= nzones; ++z) {
2975f34f2dcSJed Brown       PetscCallCGNS(cg_nbocos(cgid, 1, z, &nbc));
2985f34f2dcSJed Brown       for (bc = 1; bc <= nbc; ++bc) {
2995f34f2dcSJed Brown         PetscCallCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets));
3005f34f2dcSJed Brown         PetscCall(DMCreateLabel(*dm, bcname));
3015f34f2dcSJed Brown         PetscCall(DMGetLabel(*dm, bcname, &label));
3025f34f2dcSJed Brown         PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals));
3035f34f2dcSJed Brown         PetscCallCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *) normals));
3045f34f2dcSJed Brown         if (pointtype == CGNS_ENUMV(ElementRange)) {
3055f34f2dcSJed Brown           /* Range of cells: assuming half-open interval since the documentation sucks */
3065f34f2dcSJed Brown           for (c = points[0]; c < points[1]; ++c) {
3075f34f2dcSJed Brown             PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1));
3085f34f2dcSJed Brown           }
3095f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(ElementList)) {
3105f34f2dcSJed Brown           /* List of cells */
3115f34f2dcSJed Brown           for (c = 0; c < npoints; ++c) {
3125f34f2dcSJed Brown             PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1));
3135f34f2dcSJed Brown           }
3145f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(PointRange)) {
3155f34f2dcSJed Brown           CGNS_ENUMT(GridLocation_t) gridloc;
3165f34f2dcSJed Brown 
3175f34f2dcSJed Brown           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
3185f34f2dcSJed Brown           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
3195f34f2dcSJed Brown           PetscCallCGNS(cg_gridlocation_read(&gridloc));
3205f34f2dcSJed Brown           /* Range of points: assuming half-open interval since the documentation sucks */
3215f34f2dcSJed Brown           for (c = points[0]; c < points[1]; ++c) {
3225f34f2dcSJed Brown             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z-1], 1));
3235f34f2dcSJed Brown             else                               PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1));
3245f34f2dcSJed Brown           }
3255f34f2dcSJed Brown         } else if (pointtype == CGNS_ENUMV(PointList)) {
3265f34f2dcSJed Brown           CGNS_ENUMT(GridLocation_t) gridloc;
3275f34f2dcSJed Brown 
3285f34f2dcSJed Brown           /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */
3295f34f2dcSJed Brown           PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end"));
3305f34f2dcSJed Brown           PetscCallCGNS(cg_gridlocation_read(&gridloc));
3315f34f2dcSJed Brown           for (c = 0; c < npoints; ++c) {
3325f34f2dcSJed Brown             if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z-1], 1));
3335f34f2dcSJed Brown             else                               PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1));
3345f34f2dcSJed Brown           }
3355f34f2dcSJed Brown         } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype);
3365f34f2dcSJed Brown         PetscCall(PetscFree2(points, normals));
3375f34f2dcSJed Brown       }
3385f34f2dcSJed Brown     }
3395f34f2dcSJed Brown     PetscCall(PetscFree2(cellStart, vertStart));
3405f34f2dcSJed Brown   }
3415f34f2dcSJed Brown   PetscCall(DMGetNumLabels(*dm, &labelIdRange[1]));
3425f34f2dcSJed Brown   PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm));
3435f34f2dcSJed Brown 
3445f34f2dcSJed Brown   /* Create BC labels at all processes */
3455f34f2dcSJed Brown   for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) {
3465f34f2dcSJed Brown     char *labelName = buffer;
3475f34f2dcSJed Brown     size_t len = sizeof(buffer);
3485f34f2dcSJed Brown     const char *locName;
3495f34f2dcSJed Brown 
3505f34f2dcSJed Brown     if (rank == 0) {
3515f34f2dcSJed Brown       PetscCall(DMGetLabelByNum(*dm, labelId, &label));
3525f34f2dcSJed Brown       PetscCall(PetscObjectGetName((PetscObject)label, &locName));
3535f34f2dcSJed Brown       PetscCall(PetscStrncpy(labelName, locName, len));
3545f34f2dcSJed Brown     }
3555f34f2dcSJed Brown     PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm));
3565f34f2dcSJed Brown     PetscCallMPI(DMCreateLabel(*dm, labelName));
3575f34f2dcSJed Brown   }
3585f34f2dcSJed Brown   PetscFunctionReturn(0);
3595f34f2dcSJed Brown }
3605f34f2dcSJed Brown 
3615f34f2dcSJed Brown // Permute plex closure ordering to CGNS
362*e853fb4cSJames Wright static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) *element_type, const int **perm)
3635f34f2dcSJed Brown {
3645f34f2dcSJed Brown   // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unst_example
3655f34f2dcSJed Brown   static const int bar_2[2] = {0, 1};
3665f34f2dcSJed Brown   static const int bar_3[3] = {1, 2, 0};
3675f34f2dcSJed Brown   static const int bar_4[4] = {2, 3, 0, 1};
3685f34f2dcSJed Brown   static const int bar_5[5] = {3, 4, 0, 1, 2};
3695f34f2dcSJed Brown   static const int tri_3[3] = {0, 1, 2};
3705f34f2dcSJed Brown   static const int tri_6[6] = {3, 4, 5, 0, 1, 2};
3715f34f2dcSJed Brown   static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0};
3725f34f2dcSJed Brown   static const int quad_4[4] = {0, 1, 2, 3};
3735f34f2dcSJed Brown   static const int quad_9[9] = {
3745f34f2dcSJed Brown     5, 6, 7, 8, // vertices
3755f34f2dcSJed Brown     1, 2, 3, 4, // edges
3765f34f2dcSJed Brown     0,          // center
3775f34f2dcSJed Brown   };
3785f34f2dcSJed Brown   static const int quad_16[] = {
3795f34f2dcSJed Brown     12, 13, 14, 15,             // vertices
3805f34f2dcSJed Brown     4, 5, 6, 7, 8, 9, 10, 11,   // edges
3815f34f2dcSJed Brown     0, 1, 3, 2,                 // centers
3825f34f2dcSJed Brown   };
3835f34f2dcSJed Brown   static const int quad_25[] = {
3845f34f2dcSJed Brown     21, 22, 23, 24,                                // vertices
3855f34f2dcSJed Brown     9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges
3865f34f2dcSJed Brown     0, 1, 2, 5, 8, 7, 6, 3, 4,                     // centers
3875f34f2dcSJed Brown   };
3885f34f2dcSJed Brown   static const int tetra_4[4] = {0, 2, 1, 3};
3895f34f2dcSJed Brown   static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4};
3905f34f2dcSJed Brown   static const int tetra_20[20] = {
3915f34f2dcSJed Brown     16, 18, 17, 19,             // vertices
3925f34f2dcSJed Brown     9, 8, 7, 6, 5, 4,           // bottom edges
3935f34f2dcSJed Brown     10, 11, 14, 15, 13, 12,     // side edges
3945f34f2dcSJed Brown     0, 2, 3, 1,                 // faces
3955f34f2dcSJed Brown   };
3965f34f2dcSJed Brown   static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7};
3975f34f2dcSJed Brown   static const int hexa_27[27] = {
3985f34f2dcSJed Brown     19, 22, 21, 20, 23, 24, 25, 26, // vertices
3995f34f2dcSJed Brown     10, 9, 8, 7, // bottom edges
4005f34f2dcSJed Brown     16, 15, 18, 17, // mid edges
4015f34f2dcSJed Brown     11, 12, 13, 14, // top edges
4025f34f2dcSJed Brown     1, 3, 5, 4, 6, 2, // faces
4035f34f2dcSJed Brown     0, // center
4045f34f2dcSJed Brown   };
4055f34f2dcSJed Brown   static const int hexa_64[64] = { // 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
4065f34f2dcSJed Brown     56, 59, 58, 57, 60, 61, 62, 63, // vertices
4075f34f2dcSJed Brown     39, 38, 37, 36, 35, 34, 33, 32, // bottom edges
4085f34f2dcSJed Brown     51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24
4095f34f2dcSJed Brown     40, 41, 42, 43, 44, 45, 46, 47, // top edges
4105f34f2dcSJed Brown     8, 10, 11, 9, // z-minus face
4115f34f2dcSJed Brown     16, 17, 19, 18, // y-minus face
4125f34f2dcSJed Brown     24, 25, 27, 26, // x-plus face
4135f34f2dcSJed Brown     20, 21, 23, 22, // y-plus face
4145f34f2dcSJed Brown     30, 28, 29, 31, // x-minus face
4155f34f2dcSJed Brown     12, 13, 15, 14, // z-plus face
4165f34f2dcSJed Brown     0, 1, 3, 2, 4, 5, 7, 6, // center
4175f34f2dcSJed Brown   };
4185f34f2dcSJed Brown 
4195f34f2dcSJed Brown   PetscFunctionBegin;
4205f34f2dcSJed Brown   *element_type = CGNS_ENUMV(ElementTypeNull);
4215f34f2dcSJed Brown   *perm = NULL;
4225f34f2dcSJed Brown   switch (cell_type) {
4235f34f2dcSJed Brown   case DM_POLYTOPE_SEGMENT:
4245f34f2dcSJed Brown     switch (closure_size) {
4255f34f2dcSJed Brown     case 2:
4265f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_2);
4275f34f2dcSJed Brown       *perm = bar_2;
4285f34f2dcSJed Brown     case 3:
4295f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_3);
4305f34f2dcSJed Brown       *perm = bar_3;
4315f34f2dcSJed Brown     case 4:
4325f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_4);
4335f34f2dcSJed Brown       *perm = bar_4;
4345f34f2dcSJed Brown       break;
4355f34f2dcSJed Brown     case 5:
4365f34f2dcSJed Brown       *element_type = CGNS_ENUMV(BAR_5);
4375f34f2dcSJed Brown       *perm = bar_5;
4385f34f2dcSJed Brown       break;
4395f34f2dcSJed Brown     default:
4405f34f2dcSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
4415f34f2dcSJed Brown     }
4425f34f2dcSJed Brown     break;
4435f34f2dcSJed Brown   case DM_POLYTOPE_TRIANGLE:
4445f34f2dcSJed Brown     switch (closure_size) {
4455f34f2dcSJed Brown     case 3:
4465f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_3);
4475f34f2dcSJed Brown       *perm = tri_3;
4485f34f2dcSJed Brown       break;
4495f34f2dcSJed Brown     case 6:
4505f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_6);
4515f34f2dcSJed Brown       *perm = tri_6;
4525f34f2dcSJed Brown       break;
4535f34f2dcSJed Brown     case 10:
4545f34f2dcSJed Brown       *element_type = CGNS_ENUMV(TRI_10);
4555f34f2dcSJed Brown       *perm = tri_10;
4565f34f2dcSJed Brown       break;
4575f34f2dcSJed Brown     default:
4585f34f2dcSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
4595f34f2dcSJed Brown     }
4605f34f2dcSJed Brown     break;
4615f34f2dcSJed Brown   case DM_POLYTOPE_QUADRILATERAL:
4625f34f2dcSJed Brown     switch (closure_size) {
4635f34f2dcSJed Brown     case 4:
4645f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_4);
4655f34f2dcSJed Brown       *perm = quad_4;
4665f34f2dcSJed Brown       break;
4675f34f2dcSJed Brown     case 9:
4685f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_9);
4695f34f2dcSJed Brown       *perm = quad_9;
4705f34f2dcSJed Brown       break;
4715f34f2dcSJed Brown     case 16:
4725f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_16);
4735f34f2dcSJed Brown       *perm = quad_16;
4745f34f2dcSJed Brown       break;
4755f34f2dcSJed Brown     case 25:
4765f34f2dcSJed Brown       *element_type = CGNS_ENUMV(QUAD_25);
4775f34f2dcSJed Brown       *perm = quad_25;
4785f34f2dcSJed Brown       break;
4795f34f2dcSJed Brown     default:
4805f34f2dcSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
4815f34f2dcSJed Brown     }
4825f34f2dcSJed Brown     break;
4835f34f2dcSJed Brown   case DM_POLYTOPE_TETRAHEDRON:
4845f34f2dcSJed Brown     switch (closure_size) {
4855f34f2dcSJed Brown       case 4:
4865f34f2dcSJed Brown         *element_type = CGNS_ENUMV(TETRA_4);
4875f34f2dcSJed Brown         *perm = tetra_4;
4885f34f2dcSJed Brown         break;
4895f34f2dcSJed Brown       case 10:
4905f34f2dcSJed Brown         *element_type = CGNS_ENUMV(TETRA_10);
4915f34f2dcSJed Brown         *perm = tetra_10;
4925f34f2dcSJed Brown         break;
4935f34f2dcSJed Brown       case 20:
4945f34f2dcSJed Brown         *element_type = CGNS_ENUMV(TETRA_20);
4955f34f2dcSJed Brown         *perm = tetra_20;
4965f34f2dcSJed Brown         break;
4975f34f2dcSJed Brown     default:
4985f34f2dcSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
4995f34f2dcSJed Brown     }
5005f34f2dcSJed Brown     break;
5015f34f2dcSJed Brown   case DM_POLYTOPE_HEXAHEDRON:
5025f34f2dcSJed Brown     switch (closure_size) {
5035f34f2dcSJed Brown     case 8:
5045f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_8);
5055f34f2dcSJed Brown       *perm = hexa_8;
5065f34f2dcSJed Brown       break;
5075f34f2dcSJed Brown     case 27:
5085f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_27);
5095f34f2dcSJed Brown       *perm = hexa_27;
5105f34f2dcSJed Brown       break;
5115f34f2dcSJed Brown     case 64:
5125f34f2dcSJed Brown       *element_type = CGNS_ENUMV(HEXA_64);
5135f34f2dcSJed Brown       *perm = hexa_64;
5145f34f2dcSJed Brown       break;
5155f34f2dcSJed Brown     default:
5165f34f2dcSJed Brown       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5175f34f2dcSJed Brown     }
5185f34f2dcSJed Brown     break;
5195f34f2dcSJed Brown   default:
5205f34f2dcSJed Brown     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size);
5215f34f2dcSJed Brown   }
5225f34f2dcSJed Brown   PetscFunctionReturn(0);
5235f34f2dcSJed Brown }
5245f34f2dcSJed Brown 
5255f34f2dcSJed Brown // node_l2g must be freed
5265f34f2dcSJed Brown static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g)
5275f34f2dcSJed Brown {
5285f34f2dcSJed Brown   PetscSection local_section;
5295f34f2dcSJed Brown   PetscSF point_sf;
5305f34f2dcSJed Brown   PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes;
5315f34f2dcSJed Brown   PetscMPIInt comm_size;
5325f34f2dcSJed Brown   const PetscInt *ilocal, field = 0;
5335f34f2dcSJed Brown 
5345f34f2dcSJed Brown   PetscFunctionBegin;
5355f34f2dcSJed Brown   *num_local_nodes = -1;
5365f34f2dcSJed Brown   *num_global_nodes = -1;
5375f34f2dcSJed Brown   *nStart = -1;
5385f34f2dcSJed Brown   *nEnd =  -1;
5395f34f2dcSJed Brown   *node_l2g = NULL;
5405f34f2dcSJed Brown 
5415f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size));
5425f34f2dcSJed Brown   PetscCall(DMGetLocalSection(dm, &local_section));
5435f34f2dcSJed Brown   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
5445f34f2dcSJed Brown   PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd));
5455f34f2dcSJed Brown   PetscCall(DMGetPointSF(dm, &point_sf));
5465f34f2dcSJed Brown   if (comm_size == 1) nleaves = 0;
5475f34f2dcSJed Brown   else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL));
5485f34f2dcSJed Brown   PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp));
5495f34f2dcSJed Brown 
5505f34f2dcSJed Brown   PetscInt local_node = 0, owned_node = 0, owned_start = 0;
5515f34f2dcSJed Brown   for (PetscInt p=spStart,leaf=0; p<spEnd; p++) {
5525f34f2dcSJed Brown     PetscInt dof;
5535f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
5545f34f2dcSJed 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);
5555f34f2dcSJed Brown     local_node += dof / ncomp;
5565f34f2dcSJed Brown     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
5575f34f2dcSJed Brown       leaf++;
5585f34f2dcSJed Brown     } else {
5595f34f2dcSJed Brown       owned_node += dof / ncomp;
5605f34f2dcSJed Brown     }
5615f34f2dcSJed Brown   }
5625f34f2dcSJed Brown   PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
5635f34f2dcSJed Brown   PetscCall(PetscMalloc1(pEnd-pStart, &points));
5645f34f2dcSJed Brown   owned_node = 0;
5655f34f2dcSJed Brown   for (PetscInt p=spStart,leaf=0; p<spEnd; p++) {
5665f34f2dcSJed Brown     if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process
5675f34f2dcSJed Brown       points[p-pStart] = -1;
5685f34f2dcSJed Brown       leaf++;
5695f34f2dcSJed Brown       continue;
5705f34f2dcSJed Brown     }
5715f34f2dcSJed Brown     PetscInt dof, offset;
5725f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
5735f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
5745f34f2dcSJed 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);
5755f34f2dcSJed Brown     points[p-pStart] = owned_start + owned_node;
5765f34f2dcSJed Brown     owned_node += dof / ncomp;
5775f34f2dcSJed Brown   }
5785f34f2dcSJed Brown   if (comm_size > 1) {
5795f34f2dcSJed Brown     PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE));
5805f34f2dcSJed Brown     PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE));
5815f34f2dcSJed Brown   }
5825f34f2dcSJed Brown 
5835f34f2dcSJed Brown   // Set up global indices for each local node
5845f34f2dcSJed Brown   PetscCall(PetscMalloc1(local_node, &nodes));
5855f34f2dcSJed Brown   for (PetscInt p=spStart; p<spEnd; p++) {
5865f34f2dcSJed Brown     PetscInt dof, offset;
5875f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof));
5885f34f2dcSJed Brown     PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset));
5895f34f2dcSJed Brown     for (PetscInt n=0; n<dof/ncomp; n++) nodes[offset/ncomp+n] = points[p-pStart] + n;
5905f34f2dcSJed Brown   }
5915f34f2dcSJed Brown   PetscCall(PetscFree(points));
5925f34f2dcSJed Brown   *num_local_nodes = local_node;
5935f34f2dcSJed Brown   *nStart = owned_start;
5945f34f2dcSJed Brown   *nEnd =  owned_start + owned_node;
5955f34f2dcSJed Brown   PetscCallMPI(MPI_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
5965f34f2dcSJed Brown   *node_l2g = nodes;
5975f34f2dcSJed Brown   PetscFunctionReturn(0);
5985f34f2dcSJed Brown }
5995f34f2dcSJed Brown 
6005f34f2dcSJed Brown PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer)
6015f34f2dcSJed Brown {
6025f34f2dcSJed Brown   PetscViewer_CGNS *cgv = (PetscViewer_CGNS*)viewer->data;
6035f34f2dcSJed Brown   PetscInt topo_dim, coord_dim, num_global_elems;
6045f34f2dcSJed Brown   PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd;
6055f34f2dcSJed Brown   const PetscInt *node_l2g;
6065f34f2dcSJed Brown   Vec coord;
6075f34f2dcSJed Brown   DM cdm;
6085f34f2dcSJed Brown   PetscMPIInt size;
6095f34f2dcSJed Brown   const char *dm_name;
6105f34f2dcSJed Brown   int base, zone;
6115f34f2dcSJed Brown   cgsize_t isize[3];
6125f34f2dcSJed Brown 
6135f34f2dcSJed Brown   PetscFunctionBegin;
6145f34f2dcSJed Brown   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size));
6155f34f2dcSJed Brown   PetscCall(DMGetDimension(dm, &topo_dim));
6165f34f2dcSJed Brown   PetscCall(DMGetCoordinateDim(dm, &coord_dim));
6175f34f2dcSJed Brown   PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name));
6185f34f2dcSJed Brown   PetscCallCGNS(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base));
6195f34f2dcSJed Brown   PetscCallCGNS(cg_goto(cgv->file_num, base, NULL));
6205f34f2dcSJed Brown   PetscCallCGNS(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)));
6215f34f2dcSJed Brown 
6225f34f2dcSJed Brown   PetscCall(DMGetCoordinateDM(dm, &cdm));
6235f34f2dcSJed Brown   PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g));
6245f34f2dcSJed Brown   PetscCall(DMGetCoordinatesLocal(dm, &coord));
6255f34f2dcSJed Brown   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
6265f34f2dcSJed Brown   num_global_elems = cEnd - cStart;
6275f34f2dcSJed Brown   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm)));
6285f34f2dcSJed Brown   isize[0] = num_global_nodes;
6295f34f2dcSJed Brown   isize[1] = num_global_elems;
6305f34f2dcSJed Brown   isize[2] = 0;
6315f34f2dcSJed Brown   PetscCallCGNS(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone));
6325f34f2dcSJed Brown 
6335f34f2dcSJed Brown   {
6345f34f2dcSJed Brown     const PetscScalar *X;
6355f34f2dcSJed Brown     PetscScalar *x;
6365f34f2dcSJed Brown     int coord_ids[3];
6375f34f2dcSJed Brown 
6385f34f2dcSJed Brown     PetscCall(VecGetArrayRead(coord, &X));
6395f34f2dcSJed Brown     for (PetscInt d=0; d<coord_dim; d++) {
6405f34f2dcSJed Brown       const double exponents[] = {0, 1, 0, 0, 0};
6415f34f2dcSJed Brown       char coord_name[64];
6425f34f2dcSJed Brown       PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d));
6435f34f2dcSJed Brown       PetscCallCGNS(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]));
6445f34f2dcSJed Brown       PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL));
6455f34f2dcSJed Brown       PetscCallCGNS(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents));
6465f34f2dcSJed Brown     }
6475f34f2dcSJed Brown 
6485f34f2dcSJed Brown     DMPolytopeType cell_type;
6495f34f2dcSJed Brown     int section;
6505f34f2dcSJed Brown     cgsize_t e_owned, e_global, e_start, *conn = NULL;
6515f34f2dcSJed Brown     const int *perm;
652*e853fb4cSJames Wright     CGNS_ENUMT(ElementType_t) element_type= CGNS_ENUMV(ElementTypeNull);
6535f34f2dcSJed Brown     {
6545f34f2dcSJed Brown       PetscCall(PetscMalloc1(nEnd - nStart, &x));
6555f34f2dcSJed Brown       for (PetscInt d=0; d<coord_dim; d++) {
6565f34f2dcSJed Brown         for (PetscInt n=0; n<num_local_nodes; n++) {
6575f34f2dcSJed Brown           PetscInt gn = node_l2g[n];
6585f34f2dcSJed Brown           if (gn < nStart || nEnd <= gn) continue;
6595f34f2dcSJed Brown           x[gn - nStart] = X[n*coord_dim + d];
6605f34f2dcSJed Brown         }
6615f34f2dcSJed Brown         // CGNS nodes use 1-based indexing
6625f34f2dcSJed Brown         cgsize_t start = nStart + 1, end = nEnd;
6635f34f2dcSJed Brown         PetscCallCGNS(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x));
6645f34f2dcSJed Brown       }
6655f34f2dcSJed Brown       PetscCall(PetscFree(x));
6665f34f2dcSJed Brown       PetscCall(VecRestoreArrayRead(coord, &X));
6675f34f2dcSJed Brown     }
6685f34f2dcSJed Brown 
6695f34f2dcSJed Brown     PetscCall(DMPlexGetCellType(dm, cStart, &cell_type));
6705f34f2dcSJed Brown     for (PetscInt i=cStart,c=0; i<cEnd; i++) {
6715f34f2dcSJed Brown       PetscInt closure_dof, *closure_indices, elem_size;
6725f34f2dcSJed Brown       PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
6735f34f2dcSJed Brown       elem_size = closure_dof / coord_dim;
6745f34f2dcSJed Brown       if (!conn) PetscCall(PetscMalloc1((cEnd-cStart)*elem_size, &conn));
6755f34f2dcSJed Brown       PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm));
6765f34f2dcSJed Brown       for (PetscInt j=0; j<elem_size; j++) {
6775f34f2dcSJed Brown         conn[c++] = node_l2g[closure_indices[perm[j]*coord_dim] / coord_dim] + 1;
6785f34f2dcSJed Brown       }
6795f34f2dcSJed Brown       PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL));
6805f34f2dcSJed Brown     }
6815f34f2dcSJed Brown     e_owned = cEnd - cStart;
6825f34f2dcSJed Brown     PetscCallMPI(MPI_Allreduce(&e_owned, &e_global, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm)));
6835f34f2dcSJed 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);
6845f34f2dcSJed Brown     e_start = 0;
6855f34f2dcSJed Brown     PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm)));
6865f34f2dcSJed Brown     PetscCallCGNS(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, &section));
6875f34f2dcSJed Brown     PetscCallCGNS(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start+1, e_start + e_owned, conn));
6885f34f2dcSJed Brown     PetscCall(PetscFree(conn));
6895f34f2dcSJed Brown 
6905f34f2dcSJed Brown     cgv->base = base;
6915f34f2dcSJed Brown     cgv->zone = zone;
6925f34f2dcSJed Brown     cgv->node_l2g = node_l2g;
6935f34f2dcSJed Brown     cgv->num_local_nodes = num_local_nodes;
6945f34f2dcSJed Brown     cgv->nStart = nStart;
6955f34f2dcSJed Brown     cgv->nEnd = nEnd;
6965f34f2dcSJed Brown     if (1) {
6975f34f2dcSJed Brown       PetscMPIInt rank;
6985f34f2dcSJed Brown       int *efield;
6995f34f2dcSJed Brown       int sol, field;
7005f34f2dcSJed Brown       DMLabel label;
7015f34f2dcSJed Brown       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank));
7025f34f2dcSJed Brown       PetscCall(PetscMalloc1(e_owned, &efield));
7035f34f2dcSJed Brown       for (PetscInt i=0; i<e_owned; i++) efield[i] = rank;
7045f34f2dcSJed Brown       PetscCallCGNS(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol));
7055f34f2dcSJed Brown       PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field));
7065f34f2dcSJed Brown       cgsize_t start = e_start + 1, end = e_start + e_owned;
7075f34f2dcSJed Brown       PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
7085f34f2dcSJed Brown       PetscCall(DMGetLabel(dm, "Cell Sets", &label));
7095f34f2dcSJed Brown       if (label) {
7105f34f2dcSJed Brown         for (PetscInt c=cStart; c<cEnd; c++) {
7115f34f2dcSJed Brown           PetscInt value;
7125f34f2dcSJed Brown           PetscCall(DMLabelGetValue(label, c, &value));
7135f34f2dcSJed Brown           efield[c-cStart] = value;
7145f34f2dcSJed Brown         }
7155f34f2dcSJed Brown         PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field));
7165f34f2dcSJed Brown         PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield));
7175f34f2dcSJed Brown       }
7185f34f2dcSJed Brown       PetscCall(PetscFree(efield));
7195f34f2dcSJed Brown     }
7205f34f2dcSJed Brown   }
7215f34f2dcSJed Brown   PetscFunctionReturn(0);
7225f34f2dcSJed Brown }
7235f34f2dcSJed Brown 
7245f34f2dcSJed Brown static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) *cd)
7255f34f2dcSJed Brown {
7265f34f2dcSJed Brown   PetscFunctionBegin;
7275f34f2dcSJed Brown   switch (pd) {
7285f34f2dcSJed Brown   case PETSC_FLOAT: *cd = CGNS_ENUMV(RealSingle); break;
7295f34f2dcSJed Brown   case PETSC_DOUBLE: *cd = CGNS_ENUMV(RealDouble); break;
7305f34f2dcSJed Brown   case PETSC_COMPLEX: *cd = CGNS_ENUMV(ComplexDouble); break;
7315f34f2dcSJed Brown   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]);
7325f34f2dcSJed Brown   }
7335f34f2dcSJed Brown   PetscFunctionReturn(0);
7345f34f2dcSJed Brown }
7355f34f2dcSJed Brown 
7365f34f2dcSJed Brown PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer)
7375f34f2dcSJed Brown {
7385f34f2dcSJed Brown   PetscViewer_CGNS *cgv = (PetscViewer_CGNS*)viewer->data;
7395f34f2dcSJed Brown   DM dm;
7405f34f2dcSJed Brown   PetscSection section;
7415f34f2dcSJed Brown   PetscInt ncomp, time_step;
7425f34f2dcSJed Brown   PetscReal time, *time_slot;
7435f34f2dcSJed Brown   const PetscInt field = 0;
7445f34f2dcSJed Brown   const PetscScalar *v;
7455f34f2dcSJed Brown   char solution_name[PETSC_MAX_PATH_LEN];
7465f34f2dcSJed Brown   int sol;
7475f34f2dcSJed Brown 
7485f34f2dcSJed Brown   PetscFunctionBegin;
7495f34f2dcSJed Brown   PetscCall(VecGetDM(V, &dm));
7505f34f2dcSJed Brown   if (!cgv->node_l2g) PetscCall(DMView(dm, viewer));
7515f34f2dcSJed Brown   if (!cgv->nodal_field) PetscCall(PetscMalloc1(cgv->nEnd-cgv->nStart, &cgv->nodal_field));
7525f34f2dcSJed Brown   if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times));
7535f34f2dcSJed Brown 
7545f34f2dcSJed Brown   PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time));
7555f34f2dcSJed Brown   if (time_step < 0) {time_step = 0; time = 0.;}
7565f34f2dcSJed Brown   PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot));
7575f34f2dcSJed Brown   *time_slot = time;
7585f34f2dcSJed Brown   PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step));
7595f34f2dcSJed Brown   PetscCallCGNS(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, CGNS_ENUMV(Vertex), &sol));
7605f34f2dcSJed Brown   PetscCall(DMGetLocalSection(dm, &section));
7615f34f2dcSJed Brown   PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp));
7625f34f2dcSJed Brown   PetscCall(VecGetArrayRead(V, &v));
7635f34f2dcSJed Brown   for (PetscInt comp=0; comp < ncomp; comp++) {
7645f34f2dcSJed Brown     int cgfield;
7655f34f2dcSJed Brown     const char *comp_name;
7665f34f2dcSJed Brown     CGNS_ENUMT(DataType_t) datatype;
7675f34f2dcSJed Brown     PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name));
7685f34f2dcSJed Brown     PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype));
7695f34f2dcSJed Brown     PetscCallCGNS(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, comp_name, &cgfield));
7705f34f2dcSJed Brown     for (PetscInt n=0; n<cgv->num_local_nodes; n++) {
7715f34f2dcSJed Brown       PetscInt gn = cgv->node_l2g[n];
7725f34f2dcSJed Brown       if (gn < cgv->nStart || cgv->nEnd <= gn) continue;
7735f34f2dcSJed Brown       cgv->nodal_field[gn - cgv->nStart] = v[n*ncomp + comp];
7745f34f2dcSJed Brown     }
7755f34f2dcSJed Brown     // CGNS nodes use 1-based indexing
7765f34f2dcSJed Brown     cgsize_t start = cgv->nStart + 1, end = cgv->nEnd;
7775f34f2dcSJed Brown     PetscCallCGNS(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field));
7785f34f2dcSJed Brown   }
7795f34f2dcSJed Brown   PetscCall(VecRestoreArrayRead(V, &v));
7805f34f2dcSJed Brown   PetscFunctionReturn(0);
7815f34f2dcSJed Brown }
782