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