1*5f34f2dcSJed Brown #define PETSCDM_DLL 2*5f34f2dcSJed Brown #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3*5f34f2dcSJed Brown #include <petsc/private/viewercgnsimpl.h> 4*5f34f2dcSJed Brown 5*5f34f2dcSJed Brown #include <pcgnslib.h> 6*5f34f2dcSJed Brown #include <cgns_io.h> 7*5f34f2dcSJed Brown 8*5f34f2dcSJed Brown #if !defined(CGNS_ENUMT) 9*5f34f2dcSJed Brown #define CGNS_ENUMT(a) a 10*5f34f2dcSJed Brown #endif 11*5f34f2dcSJed Brown #if !defined(CGNS_ENUMV) 12*5f34f2dcSJed Brown #define CGNS_ENUMV(a) a 13*5f34f2dcSJed Brown #endif 14*5f34f2dcSJed Brown 15*5f34f2dcSJed Brown PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 16*5f34f2dcSJed Brown { 17*5f34f2dcSJed Brown PetscMPIInt rank; 18*5f34f2dcSJed Brown int cgid = -1; 19*5f34f2dcSJed Brown 20*5f34f2dcSJed Brown PetscFunctionBegin; 21*5f34f2dcSJed Brown PetscValidCharPointer(filename, 2); 22*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_rank(comm, &rank)); 23*5f34f2dcSJed Brown if (rank == 0) { 24*5f34f2dcSJed Brown PetscCallCGNS(cg_open(filename, CG_MODE_READ, &cgid)); 25*5f34f2dcSJed Brown PetscCheck(cgid > 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 26*5f34f2dcSJed Brown } 27*5f34f2dcSJed Brown PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 28*5f34f2dcSJed Brown if (rank == 0) PetscCallCGNS(cg_close(cgid)); 29*5f34f2dcSJed Brown PetscFunctionReturn(0); 30*5f34f2dcSJed Brown } 31*5f34f2dcSJed Brown 32*5f34f2dcSJed Brown PetscErrorCode DMPlexCreateCGNS_Internal(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 33*5f34f2dcSJed Brown { 34*5f34f2dcSJed Brown PetscMPIInt num_proc, rank; 35*5f34f2dcSJed Brown DM cdm; 36*5f34f2dcSJed Brown DMLabel label; 37*5f34f2dcSJed Brown PetscSection coordSection; 38*5f34f2dcSJed Brown Vec coordinates; 39*5f34f2dcSJed Brown PetscScalar *coords; 40*5f34f2dcSJed Brown PetscInt *cellStart, *vertStart, v; 41*5f34f2dcSJed Brown PetscInt labelIdRange[2], labelId; 42*5f34f2dcSJed Brown /* Read from file */ 43*5f34f2dcSJed Brown char basename[CGIO_MAX_NAME_LENGTH+1]; 44*5f34f2dcSJed Brown char buffer[CGIO_MAX_NAME_LENGTH+1]; 45*5f34f2dcSJed Brown int dim = 0, physDim = 0, coordDim =0, numVertices = 0, numCells = 0; 46*5f34f2dcSJed Brown int nzones = 0; 47*5f34f2dcSJed Brown 48*5f34f2dcSJed Brown PetscFunctionBegin; 49*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_rank(comm, &rank)); 50*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 51*5f34f2dcSJed Brown PetscCall(DMCreate(comm, dm)); 52*5f34f2dcSJed Brown PetscCall(DMSetType(*dm, DMPLEX)); 53*5f34f2dcSJed Brown 54*5f34f2dcSJed Brown /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 55*5f34f2dcSJed Brown if (rank == 0) { 56*5f34f2dcSJed Brown int nbases, z; 57*5f34f2dcSJed Brown 58*5f34f2dcSJed Brown PetscCallCGNS(cg_nbases(cgid, &nbases)); 59*5f34f2dcSJed Brown PetscCheck(nbases <= 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d",nbases); 60*5f34f2dcSJed Brown PetscCallCGNS(cg_base_read(cgid, 1, basename, &dim, &physDim)); 61*5f34f2dcSJed Brown PetscCallCGNS(cg_nzones(cgid, 1, &nzones)); 62*5f34f2dcSJed Brown PetscCall(PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart)); 63*5f34f2dcSJed Brown for (z = 1; z <= nzones; ++z) { 64*5f34f2dcSJed Brown cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 65*5f34f2dcSJed Brown 66*5f34f2dcSJed Brown PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes)); 67*5f34f2dcSJed Brown numVertices += sizes[0]; 68*5f34f2dcSJed Brown numCells += sizes[1]; 69*5f34f2dcSJed Brown cellStart[z] += sizes[1] + cellStart[z-1]; 70*5f34f2dcSJed Brown vertStart[z] += sizes[0] + vertStart[z-1]; 71*5f34f2dcSJed Brown } 72*5f34f2dcSJed Brown for (z = 1; z <= nzones; ++z) { 73*5f34f2dcSJed Brown vertStart[z] += numCells; 74*5f34f2dcSJed Brown } 75*5f34f2dcSJed Brown coordDim = dim; 76*5f34f2dcSJed Brown } 77*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm)); 78*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 79*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm)); 80*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm)); 81*5f34f2dcSJed Brown 82*5f34f2dcSJed Brown PetscCall(PetscObjectSetName((PetscObject) *dm, basename)); 83*5f34f2dcSJed Brown PetscCall(DMSetDimension(*dm, dim)); 84*5f34f2dcSJed Brown PetscCall(DMCreateLabel(*dm, "celltype")); 85*5f34f2dcSJed Brown PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices)); 86*5f34f2dcSJed Brown 87*5f34f2dcSJed Brown /* Read zone information */ 88*5f34f2dcSJed Brown if (rank == 0) { 89*5f34f2dcSJed Brown int z, c, c_loc; 90*5f34f2dcSJed Brown 91*5f34f2dcSJed Brown /* Read the cell set connectivity table and build mesh topology 92*5f34f2dcSJed Brown CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 93*5f34f2dcSJed Brown /* First set sizes */ 94*5f34f2dcSJed Brown for (z = 1, c = 0; z <= nzones; ++z) { 95*5f34f2dcSJed Brown CGNS_ENUMT(ZoneType_t) zonetype; 96*5f34f2dcSJed Brown int nsections; 97*5f34f2dcSJed Brown CGNS_ENUMT(ElementType_t) cellType; 98*5f34f2dcSJed Brown cgsize_t start, end; 99*5f34f2dcSJed Brown int nbndry, parentFlag; 100*5f34f2dcSJed Brown PetscInt numCorners; 101*5f34f2dcSJed Brown DMPolytopeType ctype; 102*5f34f2dcSJed Brown 103*5f34f2dcSJed Brown PetscCallCGNS(cg_zone_type(cgid, 1, z, &zonetype)); 104*5f34f2dcSJed Brown PetscCheck(zonetype != CGNS_ENUMV(Structured),PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS"); 105*5f34f2dcSJed Brown PetscCallCGNS(cg_nsections(cgid, 1, z, &nsections)); 106*5f34f2dcSJed Brown PetscCheck(nsections <= 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d",nsections); 107*5f34f2dcSJed Brown PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag)); 108*5f34f2dcSJed Brown /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 109*5f34f2dcSJed Brown if (cellType == CGNS_ENUMV(MIXED)) { 110*5f34f2dcSJed Brown cgsize_t elementDataSize, *elements; 111*5f34f2dcSJed Brown PetscInt off; 112*5f34f2dcSJed Brown 113*5f34f2dcSJed Brown PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize)); 114*5f34f2dcSJed Brown PetscCall(PetscMalloc1(elementDataSize, &elements)); 115*5f34f2dcSJed Brown PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL)); 116*5f34f2dcSJed Brown for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 117*5f34f2dcSJed Brown switch (elements[off]) { 118*5f34f2dcSJed Brown case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 119*5f34f2dcSJed Brown case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 120*5f34f2dcSJed Brown case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 121*5f34f2dcSJed Brown case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 122*5f34f2dcSJed Brown case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 123*5f34f2dcSJed Brown case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 124*5f34f2dcSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]); 125*5f34f2dcSJed Brown } 126*5f34f2dcSJed Brown PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 127*5f34f2dcSJed Brown PetscCall(DMPlexSetCellType(*dm, c, ctype)); 128*5f34f2dcSJed Brown off += numCorners+1; 129*5f34f2dcSJed Brown } 130*5f34f2dcSJed Brown PetscCall(PetscFree(elements)); 131*5f34f2dcSJed Brown } else { 132*5f34f2dcSJed Brown switch (cellType) { 133*5f34f2dcSJed Brown case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 134*5f34f2dcSJed Brown case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 135*5f34f2dcSJed Brown case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 136*5f34f2dcSJed Brown case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 137*5f34f2dcSJed Brown case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 138*5f34f2dcSJed Brown case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 139*5f34f2dcSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 140*5f34f2dcSJed Brown } 141*5f34f2dcSJed Brown for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 142*5f34f2dcSJed Brown PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 143*5f34f2dcSJed Brown PetscCall(DMPlexSetCellType(*dm, c, ctype)); 144*5f34f2dcSJed Brown } 145*5f34f2dcSJed Brown } 146*5f34f2dcSJed Brown } 147*5f34f2dcSJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 148*5f34f2dcSJed Brown PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 149*5f34f2dcSJed Brown } 150*5f34f2dcSJed Brown } 151*5f34f2dcSJed Brown 152*5f34f2dcSJed Brown PetscCall(DMSetUp(*dm)); 153*5f34f2dcSJed Brown 154*5f34f2dcSJed Brown PetscCall(DMCreateLabel(*dm, "zone")); 155*5f34f2dcSJed Brown if (rank == 0) { 156*5f34f2dcSJed Brown int z, c, c_loc, v_loc; 157*5f34f2dcSJed Brown 158*5f34f2dcSJed Brown PetscCall(DMGetLabel(*dm, "zone", &label)); 159*5f34f2dcSJed Brown for (z = 1, c = 0; z <= nzones; ++z) { 160*5f34f2dcSJed Brown CGNS_ENUMT(ElementType_t) cellType; 161*5f34f2dcSJed Brown cgsize_t elementDataSize, *elements, start, end; 162*5f34f2dcSJed Brown int nbndry, parentFlag; 163*5f34f2dcSJed Brown PetscInt *cone, numc, numCorners, maxCorners = 27; 164*5f34f2dcSJed Brown 165*5f34f2dcSJed Brown PetscCallCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag)); 166*5f34f2dcSJed Brown numc = end - start; 167*5f34f2dcSJed Brown /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 168*5f34f2dcSJed Brown PetscCallCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize)); 169*5f34f2dcSJed Brown PetscCall(PetscMalloc2(elementDataSize,&elements,maxCorners,&cone)); 170*5f34f2dcSJed Brown PetscCallCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL)); 171*5f34f2dcSJed Brown if (cellType == CGNS_ENUMV(MIXED)) { 172*5f34f2dcSJed Brown /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 173*5f34f2dcSJed Brown for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 174*5f34f2dcSJed Brown switch (elements[v]) { 175*5f34f2dcSJed Brown case CGNS_ENUMV(BAR_2): numCorners = 2; break; 176*5f34f2dcSJed Brown case CGNS_ENUMV(TRI_3): numCorners = 3; break; 177*5f34f2dcSJed Brown case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 178*5f34f2dcSJed Brown case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 179*5f34f2dcSJed Brown case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 180*5f34f2dcSJed Brown case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 181*5f34f2dcSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]); 182*5f34f2dcSJed Brown } 183*5f34f2dcSJed Brown ++v; 184*5f34f2dcSJed Brown for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 185*5f34f2dcSJed Brown cone[v_loc] = elements[v]+numCells-1; 186*5f34f2dcSJed Brown } 187*5f34f2dcSJed Brown PetscCall(DMPlexReorderCell(*dm, c, cone)); 188*5f34f2dcSJed Brown PetscCall(DMPlexSetCone(*dm, c, cone)); 189*5f34f2dcSJed Brown PetscCall(DMLabelSetValue(label, c, z)); 190*5f34f2dcSJed Brown } 191*5f34f2dcSJed Brown } else { 192*5f34f2dcSJed Brown switch (cellType) { 193*5f34f2dcSJed Brown case CGNS_ENUMV(BAR_2): numCorners = 2; break; 194*5f34f2dcSJed Brown case CGNS_ENUMV(TRI_3): numCorners = 3; break; 195*5f34f2dcSJed Brown case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 196*5f34f2dcSJed Brown case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 197*5f34f2dcSJed Brown case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 198*5f34f2dcSJed Brown case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 199*5f34f2dcSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 200*5f34f2dcSJed Brown } 201*5f34f2dcSJed Brown /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 202*5f34f2dcSJed Brown for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 203*5f34f2dcSJed Brown for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 204*5f34f2dcSJed Brown cone[v_loc] = elements[v]+numCells-1; 205*5f34f2dcSJed Brown } 206*5f34f2dcSJed Brown PetscCall(DMPlexReorderCell(*dm, c, cone)); 207*5f34f2dcSJed Brown PetscCall(DMPlexSetCone(*dm, c, cone)); 208*5f34f2dcSJed Brown PetscCall(DMLabelSetValue(label, c, z)); 209*5f34f2dcSJed Brown } 210*5f34f2dcSJed Brown } 211*5f34f2dcSJed Brown PetscCall(PetscFree2(elements,cone)); 212*5f34f2dcSJed Brown } 213*5f34f2dcSJed Brown } 214*5f34f2dcSJed Brown 215*5f34f2dcSJed Brown PetscCall(DMPlexSymmetrize(*dm)); 216*5f34f2dcSJed Brown PetscCall(DMPlexStratify(*dm)); 217*5f34f2dcSJed Brown if (interpolate) { 218*5f34f2dcSJed Brown DM idm; 219*5f34f2dcSJed Brown 220*5f34f2dcSJed Brown PetscCall(DMPlexInterpolate(*dm, &idm)); 221*5f34f2dcSJed Brown PetscCall(DMDestroy(dm)); 222*5f34f2dcSJed Brown *dm = idm; 223*5f34f2dcSJed Brown } 224*5f34f2dcSJed Brown 225*5f34f2dcSJed Brown /* Read coordinates */ 226*5f34f2dcSJed Brown PetscCall(DMSetCoordinateDim(*dm, coordDim)); 227*5f34f2dcSJed Brown PetscCall(DMGetCoordinateDM(*dm, &cdm)); 228*5f34f2dcSJed Brown PetscCall(DMGetLocalSection(cdm, &coordSection)); 229*5f34f2dcSJed Brown PetscCall(PetscSectionSetNumFields(coordSection, 1)); 230*5f34f2dcSJed Brown PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim)); 231*5f34f2dcSJed Brown PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 232*5f34f2dcSJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 233*5f34f2dcSJed Brown PetscCall(PetscSectionSetDof(coordSection, v, dim)); 234*5f34f2dcSJed Brown PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim)); 235*5f34f2dcSJed Brown } 236*5f34f2dcSJed Brown PetscCall(PetscSectionSetUp(coordSection)); 237*5f34f2dcSJed Brown 238*5f34f2dcSJed Brown PetscCall(DMCreateLocalVector(cdm, &coordinates)); 239*5f34f2dcSJed Brown PetscCall(VecGetArray(coordinates, &coords)); 240*5f34f2dcSJed Brown if (rank == 0) { 241*5f34f2dcSJed Brown PetscInt off = 0; 242*5f34f2dcSJed Brown float *x[3]; 243*5f34f2dcSJed Brown int z, d; 244*5f34f2dcSJed Brown 245*5f34f2dcSJed Brown PetscCall(PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2])); 246*5f34f2dcSJed Brown for (z = 1; z <= nzones; ++z) { 247*5f34f2dcSJed Brown CGNS_ENUMT(DataType_t) datatype; 248*5f34f2dcSJed Brown cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 249*5f34f2dcSJed Brown cgsize_t range_min[3] = {1, 1, 1}; 250*5f34f2dcSJed Brown cgsize_t range_max[3] = {1, 1, 1}; 251*5f34f2dcSJed Brown int ngrids, ncoords; 252*5f34f2dcSJed Brown 253*5f34f2dcSJed Brown PetscCallCGNS(cg_zone_read(cgid, 1, z, buffer, sizes)); 254*5f34f2dcSJed Brown range_max[0] = sizes[0]; 255*5f34f2dcSJed Brown PetscCallCGNS(cg_ngrids(cgid, 1, z, &ngrids)); 256*5f34f2dcSJed Brown PetscCheck(ngrids <= 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d",ngrids); 257*5f34f2dcSJed Brown PetscCallCGNS(cg_ncoords(cgid, 1, z, &ncoords)); 258*5f34f2dcSJed Brown PetscCheck(ncoords == coordDim,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d",ncoords); 259*5f34f2dcSJed Brown for (d = 0; d < coordDim; ++d) { 260*5f34f2dcSJed Brown PetscCallCGNS(cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer)); 261*5f34f2dcSJed Brown PetscCallCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d])); 262*5f34f2dcSJed Brown } 263*5f34f2dcSJed Brown if (coordDim >= 1) { 264*5f34f2dcSJed Brown for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v]; 265*5f34f2dcSJed Brown } 266*5f34f2dcSJed Brown if (coordDim >= 2) { 267*5f34f2dcSJed Brown for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v]; 268*5f34f2dcSJed Brown } 269*5f34f2dcSJed Brown if (coordDim >= 3) { 270*5f34f2dcSJed Brown for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v]; 271*5f34f2dcSJed Brown } 272*5f34f2dcSJed Brown off += sizes[0]; 273*5f34f2dcSJed Brown } 274*5f34f2dcSJed Brown PetscCall(PetscFree3(x[0],x[1],x[2])); 275*5f34f2dcSJed Brown } 276*5f34f2dcSJed Brown PetscCall(VecRestoreArray(coordinates, &coords)); 277*5f34f2dcSJed Brown 278*5f34f2dcSJed Brown PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 279*5f34f2dcSJed Brown PetscCall(VecSetBlockSize(coordinates, coordDim)); 280*5f34f2dcSJed Brown PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 281*5f34f2dcSJed Brown PetscCall(VecDestroy(&coordinates)); 282*5f34f2dcSJed Brown 283*5f34f2dcSJed Brown /* Read boundary conditions */ 284*5f34f2dcSJed Brown PetscCall(DMGetNumLabels(*dm, &labelIdRange[0])); 285*5f34f2dcSJed Brown if (rank == 0) { 286*5f34f2dcSJed Brown CGNS_ENUMT(BCType_t) bctype; 287*5f34f2dcSJed Brown CGNS_ENUMT(DataType_t) datatype; 288*5f34f2dcSJed Brown CGNS_ENUMT(PointSetType_t) pointtype; 289*5f34f2dcSJed Brown cgsize_t *points; 290*5f34f2dcSJed Brown PetscReal *normals; 291*5f34f2dcSJed Brown int normal[3]; 292*5f34f2dcSJed Brown char *bcname = buffer; 293*5f34f2dcSJed Brown cgsize_t npoints, nnormals; 294*5f34f2dcSJed Brown int z, nbc, bc, c, ndatasets; 295*5f34f2dcSJed Brown 296*5f34f2dcSJed Brown for (z = 1; z <= nzones; ++z) { 297*5f34f2dcSJed Brown PetscCallCGNS(cg_nbocos(cgid, 1, z, &nbc)); 298*5f34f2dcSJed Brown for (bc = 1; bc <= nbc; ++bc) { 299*5f34f2dcSJed Brown PetscCallCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets)); 300*5f34f2dcSJed Brown PetscCall(DMCreateLabel(*dm, bcname)); 301*5f34f2dcSJed Brown PetscCall(DMGetLabel(*dm, bcname, &label)); 302*5f34f2dcSJed Brown PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals)); 303*5f34f2dcSJed Brown PetscCallCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *) normals)); 304*5f34f2dcSJed Brown if (pointtype == CGNS_ENUMV(ElementRange)) { 305*5f34f2dcSJed Brown /* Range of cells: assuming half-open interval since the documentation sucks */ 306*5f34f2dcSJed Brown for (c = points[0]; c < points[1]; ++c) { 307*5f34f2dcSJed Brown PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1)); 308*5f34f2dcSJed Brown } 309*5f34f2dcSJed Brown } else if (pointtype == CGNS_ENUMV(ElementList)) { 310*5f34f2dcSJed Brown /* List of cells */ 311*5f34f2dcSJed Brown for (c = 0; c < npoints; ++c) { 312*5f34f2dcSJed Brown PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1)); 313*5f34f2dcSJed Brown } 314*5f34f2dcSJed Brown } else if (pointtype == CGNS_ENUMV(PointRange)) { 315*5f34f2dcSJed Brown CGNS_ENUMT(GridLocation_t) gridloc; 316*5f34f2dcSJed Brown 317*5f34f2dcSJed Brown /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 318*5f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 319*5f34f2dcSJed Brown PetscCallCGNS(cg_gridlocation_read(&gridloc)); 320*5f34f2dcSJed Brown /* Range of points: assuming half-open interval since the documentation sucks */ 321*5f34f2dcSJed Brown for (c = points[0]; c < points[1]; ++c) { 322*5f34f2dcSJed Brown if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z-1], 1)); 323*5f34f2dcSJed Brown else PetscCall(DMLabelSetValue(label, c - cellStart[z-1], 1)); 324*5f34f2dcSJed Brown } 325*5f34f2dcSJed Brown } else if (pointtype == CGNS_ENUMV(PointList)) { 326*5f34f2dcSJed Brown CGNS_ENUMT(GridLocation_t) gridloc; 327*5f34f2dcSJed Brown 328*5f34f2dcSJed Brown /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 329*5f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 330*5f34f2dcSJed Brown PetscCallCGNS(cg_gridlocation_read(&gridloc)); 331*5f34f2dcSJed Brown for (c = 0; c < npoints; ++c) { 332*5f34f2dcSJed Brown if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z-1], 1)); 333*5f34f2dcSJed Brown else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z-1], 1)); 334*5f34f2dcSJed Brown } 335*5f34f2dcSJed Brown } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype); 336*5f34f2dcSJed Brown PetscCall(PetscFree2(points, normals)); 337*5f34f2dcSJed Brown } 338*5f34f2dcSJed Brown } 339*5f34f2dcSJed Brown PetscCall(PetscFree2(cellStart, vertStart)); 340*5f34f2dcSJed Brown } 341*5f34f2dcSJed Brown PetscCall(DMGetNumLabels(*dm, &labelIdRange[1])); 342*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm)); 343*5f34f2dcSJed Brown 344*5f34f2dcSJed Brown /* Create BC labels at all processes */ 345*5f34f2dcSJed Brown for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 346*5f34f2dcSJed Brown char *labelName = buffer; 347*5f34f2dcSJed Brown size_t len = sizeof(buffer); 348*5f34f2dcSJed Brown const char *locName; 349*5f34f2dcSJed Brown 350*5f34f2dcSJed Brown if (rank == 0) { 351*5f34f2dcSJed Brown PetscCall(DMGetLabelByNum(*dm, labelId, &label)); 352*5f34f2dcSJed Brown PetscCall(PetscObjectGetName((PetscObject)label, &locName)); 353*5f34f2dcSJed Brown PetscCall(PetscStrncpy(labelName, locName, len)); 354*5f34f2dcSJed Brown } 355*5f34f2dcSJed Brown PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm)); 356*5f34f2dcSJed Brown PetscCallMPI(DMCreateLabel(*dm, labelName)); 357*5f34f2dcSJed Brown } 358*5f34f2dcSJed Brown PetscFunctionReturn(0); 359*5f34f2dcSJed Brown } 360*5f34f2dcSJed Brown 361*5f34f2dcSJed Brown // Permute plex closure ordering to CGNS 362*5f34f2dcSJed Brown static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, ElementType_t *element_type, const int **perm) 363*5f34f2dcSJed Brown { 364*5f34f2dcSJed Brown // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unst_example 365*5f34f2dcSJed Brown static const int bar_2[2] = {0, 1}; 366*5f34f2dcSJed Brown static const int bar_3[3] = {1, 2, 0}; 367*5f34f2dcSJed Brown static const int bar_4[4] = {2, 3, 0, 1}; 368*5f34f2dcSJed Brown static const int bar_5[5] = {3, 4, 0, 1, 2}; 369*5f34f2dcSJed Brown static const int tri_3[3] = {0, 1, 2}; 370*5f34f2dcSJed Brown static const int tri_6[6] = {3, 4, 5, 0, 1, 2}; 371*5f34f2dcSJed Brown static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0}; 372*5f34f2dcSJed Brown static const int quad_4[4] = {0, 1, 2, 3}; 373*5f34f2dcSJed Brown static const int quad_9[9] = { 374*5f34f2dcSJed Brown 5, 6, 7, 8, // vertices 375*5f34f2dcSJed Brown 1, 2, 3, 4, // edges 376*5f34f2dcSJed Brown 0, // center 377*5f34f2dcSJed Brown }; 378*5f34f2dcSJed Brown static const int quad_16[] = { 379*5f34f2dcSJed Brown 12, 13, 14, 15, // vertices 380*5f34f2dcSJed Brown 4, 5, 6, 7, 8, 9, 10, 11, // edges 381*5f34f2dcSJed Brown 0, 1, 3, 2, // centers 382*5f34f2dcSJed Brown }; 383*5f34f2dcSJed Brown static const int quad_25[] = { 384*5f34f2dcSJed Brown 21, 22, 23, 24, // vertices 385*5f34f2dcSJed Brown 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges 386*5f34f2dcSJed Brown 0, 1, 2, 5, 8, 7, 6, 3, 4, // centers 387*5f34f2dcSJed Brown }; 388*5f34f2dcSJed Brown static const int tetra_4[4] = {0, 2, 1, 3}; 389*5f34f2dcSJed Brown static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4}; 390*5f34f2dcSJed Brown static const int tetra_20[20] = { 391*5f34f2dcSJed Brown 16, 18, 17, 19, // vertices 392*5f34f2dcSJed Brown 9, 8, 7, 6, 5, 4, // bottom edges 393*5f34f2dcSJed Brown 10, 11, 14, 15, 13, 12, // side edges 394*5f34f2dcSJed Brown 0, 2, 3, 1, // faces 395*5f34f2dcSJed Brown }; 396*5f34f2dcSJed Brown static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 397*5f34f2dcSJed Brown static const int hexa_27[27] = { 398*5f34f2dcSJed Brown 19, 22, 21, 20, 23, 24, 25, 26, // vertices 399*5f34f2dcSJed Brown 10, 9, 8, 7, // bottom edges 400*5f34f2dcSJed Brown 16, 15, 18, 17, // mid edges 401*5f34f2dcSJed Brown 11, 12, 13, 14, // top edges 402*5f34f2dcSJed Brown 1, 3, 5, 4, 6, 2, // faces 403*5f34f2dcSJed Brown 0, // center 404*5f34f2dcSJed Brown }; 405*5f34f2dcSJed 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 406*5f34f2dcSJed Brown 56, 59, 58, 57, 60, 61, 62, 63, // vertices 407*5f34f2dcSJed Brown 39, 38, 37, 36, 35, 34, 33, 32, // bottom edges 408*5f34f2dcSJed Brown 51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24 409*5f34f2dcSJed Brown 40, 41, 42, 43, 44, 45, 46, 47, // top edges 410*5f34f2dcSJed Brown 8, 10, 11, 9, // z-minus face 411*5f34f2dcSJed Brown 16, 17, 19, 18, // y-minus face 412*5f34f2dcSJed Brown 24, 25, 27, 26, // x-plus face 413*5f34f2dcSJed Brown 20, 21, 23, 22, // y-plus face 414*5f34f2dcSJed Brown 30, 28, 29, 31, // x-minus face 415*5f34f2dcSJed Brown 12, 13, 15, 14, // z-plus face 416*5f34f2dcSJed Brown 0, 1, 3, 2, 4, 5, 7, 6, // center 417*5f34f2dcSJed Brown }; 418*5f34f2dcSJed Brown 419*5f34f2dcSJed Brown PetscFunctionBegin; 420*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(ElementTypeNull); 421*5f34f2dcSJed Brown *perm = NULL; 422*5f34f2dcSJed Brown switch (cell_type) { 423*5f34f2dcSJed Brown case DM_POLYTOPE_SEGMENT: 424*5f34f2dcSJed Brown switch (closure_size) { 425*5f34f2dcSJed Brown case 2: 426*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(BAR_2); 427*5f34f2dcSJed Brown *perm = bar_2; 428*5f34f2dcSJed Brown case 3: 429*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(BAR_3); 430*5f34f2dcSJed Brown *perm = bar_3; 431*5f34f2dcSJed Brown case 4: 432*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(BAR_4); 433*5f34f2dcSJed Brown *perm = bar_4; 434*5f34f2dcSJed Brown break; 435*5f34f2dcSJed Brown case 5: 436*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(BAR_5); 437*5f34f2dcSJed Brown *perm = bar_5; 438*5f34f2dcSJed Brown break; 439*5f34f2dcSJed Brown default: 440*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 441*5f34f2dcSJed Brown } 442*5f34f2dcSJed Brown break; 443*5f34f2dcSJed Brown case DM_POLYTOPE_TRIANGLE: 444*5f34f2dcSJed Brown switch (closure_size) { 445*5f34f2dcSJed Brown case 3: 446*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TRI_3); 447*5f34f2dcSJed Brown *perm = tri_3; 448*5f34f2dcSJed Brown break; 449*5f34f2dcSJed Brown case 6: 450*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TRI_6); 451*5f34f2dcSJed Brown *perm = tri_6; 452*5f34f2dcSJed Brown break; 453*5f34f2dcSJed Brown case 10: 454*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TRI_10); 455*5f34f2dcSJed Brown *perm = tri_10; 456*5f34f2dcSJed Brown break; 457*5f34f2dcSJed Brown default: 458*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 459*5f34f2dcSJed Brown } 460*5f34f2dcSJed Brown break; 461*5f34f2dcSJed Brown case DM_POLYTOPE_QUADRILATERAL: 462*5f34f2dcSJed Brown switch (closure_size) { 463*5f34f2dcSJed Brown case 4: 464*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(QUAD_4); 465*5f34f2dcSJed Brown *perm = quad_4; 466*5f34f2dcSJed Brown break; 467*5f34f2dcSJed Brown case 9: 468*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(QUAD_9); 469*5f34f2dcSJed Brown *perm = quad_9; 470*5f34f2dcSJed Brown break; 471*5f34f2dcSJed Brown case 16: 472*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(QUAD_16); 473*5f34f2dcSJed Brown *perm = quad_16; 474*5f34f2dcSJed Brown break; 475*5f34f2dcSJed Brown case 25: 476*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(QUAD_25); 477*5f34f2dcSJed Brown *perm = quad_25; 478*5f34f2dcSJed Brown break; 479*5f34f2dcSJed Brown default: 480*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 481*5f34f2dcSJed Brown } 482*5f34f2dcSJed Brown break; 483*5f34f2dcSJed Brown case DM_POLYTOPE_TETRAHEDRON: 484*5f34f2dcSJed Brown switch (closure_size) { 485*5f34f2dcSJed Brown case 4: 486*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TETRA_4); 487*5f34f2dcSJed Brown *perm = tetra_4; 488*5f34f2dcSJed Brown break; 489*5f34f2dcSJed Brown case 10: 490*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TETRA_10); 491*5f34f2dcSJed Brown *perm = tetra_10; 492*5f34f2dcSJed Brown break; 493*5f34f2dcSJed Brown case 20: 494*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(TETRA_20); 495*5f34f2dcSJed Brown *perm = tetra_20; 496*5f34f2dcSJed Brown break; 497*5f34f2dcSJed Brown default: 498*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 499*5f34f2dcSJed Brown } 500*5f34f2dcSJed Brown break; 501*5f34f2dcSJed Brown case DM_POLYTOPE_HEXAHEDRON: 502*5f34f2dcSJed Brown switch (closure_size) { 503*5f34f2dcSJed Brown case 8: 504*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(HEXA_8); 505*5f34f2dcSJed Brown *perm = hexa_8; 506*5f34f2dcSJed Brown break; 507*5f34f2dcSJed Brown case 27: 508*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(HEXA_27); 509*5f34f2dcSJed Brown *perm = hexa_27; 510*5f34f2dcSJed Brown break; 511*5f34f2dcSJed Brown case 64: 512*5f34f2dcSJed Brown *element_type = CGNS_ENUMV(HEXA_64); 513*5f34f2dcSJed Brown *perm = hexa_64; 514*5f34f2dcSJed Brown break; 515*5f34f2dcSJed Brown default: 516*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 517*5f34f2dcSJed Brown } 518*5f34f2dcSJed Brown break; 519*5f34f2dcSJed Brown default: 520*5f34f2dcSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 521*5f34f2dcSJed Brown } 522*5f34f2dcSJed Brown PetscFunctionReturn(0); 523*5f34f2dcSJed Brown } 524*5f34f2dcSJed Brown 525*5f34f2dcSJed Brown // node_l2g must be freed 526*5f34f2dcSJed Brown static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g) 527*5f34f2dcSJed Brown { 528*5f34f2dcSJed Brown PetscSection local_section; 529*5f34f2dcSJed Brown PetscSF point_sf; 530*5f34f2dcSJed Brown PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes; 531*5f34f2dcSJed Brown PetscMPIInt comm_size; 532*5f34f2dcSJed Brown const PetscInt *ilocal, field = 0; 533*5f34f2dcSJed Brown 534*5f34f2dcSJed Brown PetscFunctionBegin; 535*5f34f2dcSJed Brown *num_local_nodes = -1; 536*5f34f2dcSJed Brown *num_global_nodes = -1; 537*5f34f2dcSJed Brown *nStart = -1; 538*5f34f2dcSJed Brown *nEnd = -1; 539*5f34f2dcSJed Brown *node_l2g = NULL; 540*5f34f2dcSJed Brown 541*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size)); 542*5f34f2dcSJed Brown PetscCall(DMGetLocalSection(dm, &local_section)); 543*5f34f2dcSJed Brown PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 544*5f34f2dcSJed Brown PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd)); 545*5f34f2dcSJed Brown PetscCall(DMGetPointSF(dm, &point_sf)); 546*5f34f2dcSJed Brown if (comm_size == 1) nleaves = 0; 547*5f34f2dcSJed Brown else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL)); 548*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp)); 549*5f34f2dcSJed Brown 550*5f34f2dcSJed Brown PetscInt local_node = 0, owned_node = 0, owned_start = 0; 551*5f34f2dcSJed Brown for (PetscInt p=spStart,leaf=0; p<spEnd; p++) { 552*5f34f2dcSJed Brown PetscInt dof; 553*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 554*5f34f2dcSJed 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); 555*5f34f2dcSJed Brown local_node += dof / ncomp; 556*5f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 557*5f34f2dcSJed Brown leaf++; 558*5f34f2dcSJed Brown } else { 559*5f34f2dcSJed Brown owned_node += dof / ncomp; 560*5f34f2dcSJed Brown } 561*5f34f2dcSJed Brown } 562*5f34f2dcSJed Brown PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 563*5f34f2dcSJed Brown PetscCall(PetscMalloc1(pEnd-pStart, &points)); 564*5f34f2dcSJed Brown owned_node = 0; 565*5f34f2dcSJed Brown for (PetscInt p=spStart,leaf=0; p<spEnd; p++) { 566*5f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 567*5f34f2dcSJed Brown points[p-pStart] = -1; 568*5f34f2dcSJed Brown leaf++; 569*5f34f2dcSJed Brown continue; 570*5f34f2dcSJed Brown } 571*5f34f2dcSJed Brown PetscInt dof, offset; 572*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 573*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 574*5f34f2dcSJed 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); 575*5f34f2dcSJed Brown points[p-pStart] = owned_start + owned_node; 576*5f34f2dcSJed Brown owned_node += dof / ncomp; 577*5f34f2dcSJed Brown } 578*5f34f2dcSJed Brown if (comm_size > 1) { 579*5f34f2dcSJed Brown PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 580*5f34f2dcSJed Brown PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 581*5f34f2dcSJed Brown } 582*5f34f2dcSJed Brown 583*5f34f2dcSJed Brown // Set up global indices for each local node 584*5f34f2dcSJed Brown PetscCall(PetscMalloc1(local_node, &nodes)); 585*5f34f2dcSJed Brown for (PetscInt p=spStart; p<spEnd; p++) { 586*5f34f2dcSJed Brown PetscInt dof, offset; 587*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 588*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 589*5f34f2dcSJed Brown for (PetscInt n=0; n<dof/ncomp; n++) nodes[offset/ncomp+n] = points[p-pStart] + n; 590*5f34f2dcSJed Brown } 591*5f34f2dcSJed Brown PetscCall(PetscFree(points)); 592*5f34f2dcSJed Brown *num_local_nodes = local_node; 593*5f34f2dcSJed Brown *nStart = owned_start; 594*5f34f2dcSJed Brown *nEnd = owned_start + owned_node; 595*5f34f2dcSJed Brown PetscCallMPI(MPI_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 596*5f34f2dcSJed Brown *node_l2g = nodes; 597*5f34f2dcSJed Brown PetscFunctionReturn(0); 598*5f34f2dcSJed Brown } 599*5f34f2dcSJed Brown 600*5f34f2dcSJed Brown PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer) 601*5f34f2dcSJed Brown { 602*5f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS*)viewer->data; 603*5f34f2dcSJed Brown PetscInt topo_dim, coord_dim, num_global_elems; 604*5f34f2dcSJed Brown PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd; 605*5f34f2dcSJed Brown const PetscInt *node_l2g; 606*5f34f2dcSJed Brown Vec coord; 607*5f34f2dcSJed Brown DM cdm; 608*5f34f2dcSJed Brown PetscMPIInt size; 609*5f34f2dcSJed Brown const char *dm_name; 610*5f34f2dcSJed Brown int base, zone; 611*5f34f2dcSJed Brown cgsize_t isize[3]; 612*5f34f2dcSJed Brown 613*5f34f2dcSJed Brown PetscFunctionBegin; 614*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 615*5f34f2dcSJed Brown PetscCall(DMGetDimension(dm, &topo_dim)); 616*5f34f2dcSJed Brown PetscCall(DMGetCoordinateDim(dm, &coord_dim)); 617*5f34f2dcSJed Brown PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name)); 618*5f34f2dcSJed Brown PetscCallCGNS(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base)); 619*5f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, NULL)); 620*5f34f2dcSJed Brown PetscCallCGNS(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional))); 621*5f34f2dcSJed Brown 622*5f34f2dcSJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 623*5f34f2dcSJed Brown PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g)); 624*5f34f2dcSJed Brown PetscCall(DMGetCoordinatesLocal(dm, &coord)); 625*5f34f2dcSJed Brown PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 626*5f34f2dcSJed Brown num_global_elems = cEnd - cStart; 627*5f34f2dcSJed Brown PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 628*5f34f2dcSJed Brown isize[0] = num_global_nodes; 629*5f34f2dcSJed Brown isize[1] = num_global_elems; 630*5f34f2dcSJed Brown isize[2] = 0; 631*5f34f2dcSJed Brown PetscCallCGNS(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone)); 632*5f34f2dcSJed Brown 633*5f34f2dcSJed Brown { 634*5f34f2dcSJed Brown const PetscScalar *X; 635*5f34f2dcSJed Brown PetscScalar *x; 636*5f34f2dcSJed Brown int coord_ids[3]; 637*5f34f2dcSJed Brown 638*5f34f2dcSJed Brown PetscCall(VecGetArrayRead(coord, &X)); 639*5f34f2dcSJed Brown for (PetscInt d=0; d<coord_dim; d++) { 640*5f34f2dcSJed Brown const double exponents[] = {0, 1, 0, 0, 0}; 641*5f34f2dcSJed Brown char coord_name[64]; 642*5f34f2dcSJed Brown PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d)); 643*5f34f2dcSJed Brown PetscCallCGNS(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d])); 644*5f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL)); 645*5f34f2dcSJed Brown PetscCallCGNS(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents)); 646*5f34f2dcSJed Brown } 647*5f34f2dcSJed Brown 648*5f34f2dcSJed Brown DMPolytopeType cell_type; 649*5f34f2dcSJed Brown int section; 650*5f34f2dcSJed Brown cgsize_t e_owned, e_global, e_start, *conn = NULL; 651*5f34f2dcSJed Brown const int *perm; 652*5f34f2dcSJed Brown ElementType_t element_type= CGNS_ENUMV(ElementTypeNull); 653*5f34f2dcSJed Brown { 654*5f34f2dcSJed Brown PetscCall(PetscMalloc1(nEnd - nStart, &x)); 655*5f34f2dcSJed Brown for (PetscInt d=0; d<coord_dim; d++) { 656*5f34f2dcSJed Brown for (PetscInt n=0; n<num_local_nodes; n++) { 657*5f34f2dcSJed Brown PetscInt gn = node_l2g[n]; 658*5f34f2dcSJed Brown if (gn < nStart || nEnd <= gn) continue; 659*5f34f2dcSJed Brown x[gn - nStart] = X[n*coord_dim + d]; 660*5f34f2dcSJed Brown } 661*5f34f2dcSJed Brown // CGNS nodes use 1-based indexing 662*5f34f2dcSJed Brown cgsize_t start = nStart + 1, end = nEnd; 663*5f34f2dcSJed Brown PetscCallCGNS(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x)); 664*5f34f2dcSJed Brown } 665*5f34f2dcSJed Brown PetscCall(PetscFree(x)); 666*5f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(coord, &X)); 667*5f34f2dcSJed Brown } 668*5f34f2dcSJed Brown 669*5f34f2dcSJed Brown PetscCall(DMPlexGetCellType(dm, cStart, &cell_type)); 670*5f34f2dcSJed Brown for (PetscInt i=cStart,c=0; i<cEnd; i++) { 671*5f34f2dcSJed Brown PetscInt closure_dof, *closure_indices, elem_size; 672*5f34f2dcSJed Brown PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 673*5f34f2dcSJed Brown elem_size = closure_dof / coord_dim; 674*5f34f2dcSJed Brown if (!conn) PetscCall(PetscMalloc1((cEnd-cStart)*elem_size, &conn)); 675*5f34f2dcSJed Brown PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm)); 676*5f34f2dcSJed Brown for (PetscInt j=0; j<elem_size; j++) { 677*5f34f2dcSJed Brown conn[c++] = node_l2g[closure_indices[perm[j]*coord_dim] / coord_dim] + 1; 678*5f34f2dcSJed Brown } 679*5f34f2dcSJed Brown PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 680*5f34f2dcSJed Brown } 681*5f34f2dcSJed Brown e_owned = cEnd - cStart; 682*5f34f2dcSJed Brown PetscCallMPI(MPI_Allreduce(&e_owned, &e_global, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm))); 683*5f34f2dcSJed 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); 684*5f34f2dcSJed Brown e_start = 0; 685*5f34f2dcSJed Brown PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_INT64, MPI_SUM, PetscObjectComm((PetscObject)dm))); 686*5f34f2dcSJed Brown PetscCallCGNS(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, §ion)); 687*5f34f2dcSJed Brown PetscCallCGNS(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start+1, e_start + e_owned, conn)); 688*5f34f2dcSJed Brown PetscCall(PetscFree(conn)); 689*5f34f2dcSJed Brown 690*5f34f2dcSJed Brown cgv->base = base; 691*5f34f2dcSJed Brown cgv->zone = zone; 692*5f34f2dcSJed Brown cgv->node_l2g = node_l2g; 693*5f34f2dcSJed Brown cgv->num_local_nodes = num_local_nodes; 694*5f34f2dcSJed Brown cgv->nStart = nStart; 695*5f34f2dcSJed Brown cgv->nEnd = nEnd; 696*5f34f2dcSJed Brown if (1) { 697*5f34f2dcSJed Brown PetscMPIInt rank; 698*5f34f2dcSJed Brown int *efield; 699*5f34f2dcSJed Brown int sol, field; 700*5f34f2dcSJed Brown DMLabel label; 701*5f34f2dcSJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 702*5f34f2dcSJed Brown PetscCall(PetscMalloc1(e_owned, &efield)); 703*5f34f2dcSJed Brown for (PetscInt i=0; i<e_owned; i++) efield[i] = rank; 704*5f34f2dcSJed Brown PetscCallCGNS(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol)); 705*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field)); 706*5f34f2dcSJed Brown cgsize_t start = e_start + 1, end = e_start + e_owned; 707*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield)); 708*5f34f2dcSJed Brown PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 709*5f34f2dcSJed Brown if (label) { 710*5f34f2dcSJed Brown for (PetscInt c=cStart; c<cEnd; c++) { 711*5f34f2dcSJed Brown PetscInt value; 712*5f34f2dcSJed Brown PetscCall(DMLabelGetValue(label, c, &value)); 713*5f34f2dcSJed Brown efield[c-cStart] = value; 714*5f34f2dcSJed Brown } 715*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field)); 716*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield)); 717*5f34f2dcSJed Brown } 718*5f34f2dcSJed Brown PetscCall(PetscFree(efield)); 719*5f34f2dcSJed Brown } 720*5f34f2dcSJed Brown } 721*5f34f2dcSJed Brown PetscFunctionReturn(0); 722*5f34f2dcSJed Brown } 723*5f34f2dcSJed Brown 724*5f34f2dcSJed Brown static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) *cd) 725*5f34f2dcSJed Brown { 726*5f34f2dcSJed Brown PetscFunctionBegin; 727*5f34f2dcSJed Brown switch (pd) { 728*5f34f2dcSJed Brown case PETSC_FLOAT: *cd = CGNS_ENUMV(RealSingle); break; 729*5f34f2dcSJed Brown case PETSC_DOUBLE: *cd = CGNS_ENUMV(RealDouble); break; 730*5f34f2dcSJed Brown case PETSC_COMPLEX: *cd = CGNS_ENUMV(ComplexDouble); break; 731*5f34f2dcSJed Brown default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]); 732*5f34f2dcSJed Brown } 733*5f34f2dcSJed Brown PetscFunctionReturn(0); 734*5f34f2dcSJed Brown } 735*5f34f2dcSJed Brown 736*5f34f2dcSJed Brown PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer) 737*5f34f2dcSJed Brown { 738*5f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS*)viewer->data; 739*5f34f2dcSJed Brown DM dm; 740*5f34f2dcSJed Brown PetscSection section; 741*5f34f2dcSJed Brown PetscInt ncomp, time_step; 742*5f34f2dcSJed Brown PetscReal time, *time_slot; 743*5f34f2dcSJed Brown const PetscInt field = 0; 744*5f34f2dcSJed Brown const PetscScalar *v; 745*5f34f2dcSJed Brown char solution_name[PETSC_MAX_PATH_LEN]; 746*5f34f2dcSJed Brown int sol; 747*5f34f2dcSJed Brown 748*5f34f2dcSJed Brown PetscFunctionBegin; 749*5f34f2dcSJed Brown PetscCall(VecGetDM(V, &dm)); 750*5f34f2dcSJed Brown if (!cgv->node_l2g) PetscCall(DMView(dm, viewer)); 751*5f34f2dcSJed Brown if (!cgv->nodal_field) PetscCall(PetscMalloc1(cgv->nEnd-cgv->nStart, &cgv->nodal_field)); 752*5f34f2dcSJed Brown if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times)); 753*5f34f2dcSJed Brown 754*5f34f2dcSJed Brown PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time)); 755*5f34f2dcSJed Brown if (time_step < 0) {time_step = 0; time = 0.;} 756*5f34f2dcSJed Brown PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot)); 757*5f34f2dcSJed Brown *time_slot = time; 758*5f34f2dcSJed Brown PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step)); 759*5f34f2dcSJed Brown PetscCallCGNS(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, CGNS_ENUMV(Vertex), &sol)); 760*5f34f2dcSJed Brown PetscCall(DMGetLocalSection(dm, §ion)); 761*5f34f2dcSJed Brown PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp)); 762*5f34f2dcSJed Brown PetscCall(VecGetArrayRead(V, &v)); 763*5f34f2dcSJed Brown for (PetscInt comp=0; comp < ncomp; comp++) { 764*5f34f2dcSJed Brown int cgfield; 765*5f34f2dcSJed Brown const char *comp_name; 766*5f34f2dcSJed Brown CGNS_ENUMT(DataType_t) datatype; 767*5f34f2dcSJed Brown PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name)); 768*5f34f2dcSJed Brown PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype)); 769*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, comp_name, &cgfield)); 770*5f34f2dcSJed Brown for (PetscInt n=0; n<cgv->num_local_nodes; n++) { 771*5f34f2dcSJed Brown PetscInt gn = cgv->node_l2g[n]; 772*5f34f2dcSJed Brown if (gn < cgv->nStart || cgv->nEnd <= gn) continue; 773*5f34f2dcSJed Brown cgv->nodal_field[gn - cgv->nStart] = v[n*ncomp + comp]; 774*5f34f2dcSJed Brown } 775*5f34f2dcSJed Brown // CGNS nodes use 1-based indexing 776*5f34f2dcSJed Brown cgsize_t start = cgv->nStart + 1, end = cgv->nEnd; 777*5f34f2dcSJed Brown PetscCallCGNS(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field)); 778*5f34f2dcSJed Brown } 779*5f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(V, &v)); 780*5f34f2dcSJed Brown PetscFunctionReturn(0); 781*5f34f2dcSJed Brown } 782