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