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 // Permute plex closure ordering to CGNS 15 static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm) 16 { 17 CGNS_ENUMT(ElementType_t) element_type_tmp; 18 19 // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 20 static const int bar_2[2] = {0, 1}; 21 static const int bar_3[3] = {1, 2, 0}; 22 static const int bar_4[4] = {2, 3, 0, 1}; 23 static const int bar_5[5] = {3, 4, 0, 1, 2}; 24 static const int tri_3[3] = {0, 1, 2}; 25 static const int tri_6[6] = {3, 4, 5, 0, 1, 2}; 26 static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0}; 27 static const int quad_4[4] = {0, 1, 2, 3}; 28 static const int quad_9[9] = { 29 5, 6, 7, 8, // vertices 30 1, 2, 3, 4, // edges 31 0, // center 32 }; 33 static const int quad_16[] = { 34 12, 13, 14, 15, // vertices 35 4, 5, 6, 7, 8, 9, 10, 11, // edges 36 0, 1, 3, 2, // centers 37 }; 38 static const int quad_25[] = { 39 21, 22, 23, 24, // vertices 40 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges 41 0, 1, 2, 5, 8, 7, 6, 3, 4, // centers 42 }; 43 static const int tetra_4[4] = {0, 2, 1, 3}; 44 static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4}; 45 static const int tetra_20[20] = { 46 16, 18, 17, 19, // vertices 47 9, 8, 7, 6, 5, 4, // bottom edges 48 10, 11, 14, 15, 13, 12, // side edges 49 0, 2, 3, 1, // faces 50 }; 51 static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 52 static const int hexa_27[27] = { 53 19, 22, 21, 20, 23, 24, 25, 26, // vertices 54 10, 9, 8, 7, // bottom edges 55 16, 15, 18, 17, // mid edges 56 11, 12, 13, 14, // top edges 57 1, 3, 5, 4, 6, 2, // faces 58 0, // center 59 }; 60 static const int hexa_64[64] = { 61 // 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 62 56, 59, 58, 57, 60, 61, 62, 63, // vertices 63 39, 38, 37, 36, 35, 34, 33, 32, // bottom edges 64 51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24 65 40, 41, 42, 43, 44, 45, 46, 47, // top edges 66 8, 10, 11, 9, // z-minus face 67 16, 17, 19, 18, // y-minus face 68 24, 25, 27, 26, // x-plus face 69 20, 21, 23, 22, // y-plus face 70 30, 28, 29, 31, // x-minus face 71 12, 13, 15, 14, // z-plus face 72 0, 1, 3, 2, 4, 5, 7, 6, // center 73 }; 74 75 PetscFunctionBegin; 76 element_type_tmp = CGNS_ENUMV(ElementTypeNull); 77 *perm = NULL; 78 switch (cell_type) { 79 case DM_POLYTOPE_SEGMENT: 80 switch (closure_size) { 81 case 2: 82 element_type_tmp = CGNS_ENUMV(BAR_2); 83 *perm = bar_2; 84 case 3: 85 element_type_tmp = CGNS_ENUMV(BAR_3); 86 *perm = bar_3; 87 case 4: 88 element_type_tmp = CGNS_ENUMV(BAR_4); 89 *perm = bar_4; 90 break; 91 case 5: 92 element_type_tmp = CGNS_ENUMV(BAR_5); 93 *perm = bar_5; 94 break; 95 default: 96 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 97 } 98 break; 99 case DM_POLYTOPE_TRIANGLE: 100 switch (closure_size) { 101 case 3: 102 element_type_tmp = CGNS_ENUMV(TRI_3); 103 *perm = tri_3; 104 break; 105 case 6: 106 element_type_tmp = CGNS_ENUMV(TRI_6); 107 *perm = tri_6; 108 break; 109 case 10: 110 element_type_tmp = CGNS_ENUMV(TRI_10); 111 *perm = tri_10; 112 break; 113 default: 114 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 115 } 116 break; 117 case DM_POLYTOPE_QUADRILATERAL: 118 switch (closure_size) { 119 case 4: 120 element_type_tmp = CGNS_ENUMV(QUAD_4); 121 *perm = quad_4; 122 break; 123 case 9: 124 element_type_tmp = CGNS_ENUMV(QUAD_9); 125 *perm = quad_9; 126 break; 127 case 16: 128 element_type_tmp = CGNS_ENUMV(QUAD_16); 129 *perm = quad_16; 130 break; 131 case 25: 132 element_type_tmp = CGNS_ENUMV(QUAD_25); 133 *perm = quad_25; 134 break; 135 default: 136 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 137 } 138 break; 139 case DM_POLYTOPE_TETRAHEDRON: 140 switch (closure_size) { 141 case 4: 142 element_type_tmp = CGNS_ENUMV(TETRA_4); 143 *perm = tetra_4; 144 break; 145 case 10: 146 element_type_tmp = CGNS_ENUMV(TETRA_10); 147 *perm = tetra_10; 148 break; 149 case 20: 150 element_type_tmp = CGNS_ENUMV(TETRA_20); 151 *perm = tetra_20; 152 break; 153 default: 154 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 155 } 156 break; 157 case DM_POLYTOPE_HEXAHEDRON: 158 switch (closure_size) { 159 case 8: 160 element_type_tmp = CGNS_ENUMV(HEXA_8); 161 *perm = hexa_8; 162 break; 163 case 27: 164 element_type_tmp = CGNS_ENUMV(HEXA_27); 165 *perm = hexa_27; 166 break; 167 case 64: 168 element_type_tmp = CGNS_ENUMV(HEXA_64); 169 *perm = hexa_64; 170 break; 171 default: 172 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 173 } 174 break; 175 default: 176 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 177 } 178 if (element_type) *element_type = element_type_tmp; 179 PetscFunctionReturn(PETSC_SUCCESS); 180 } 181 182 /* 183 Input Parameters: 184 + cellType - The CGNS-defined element type 185 186 Output Parameters: 187 + dmcelltype - The equivalent DMPolytopeType for the cellType 188 . numCorners - Number of corners of the polytope 189 - dim - The topological dimension of the polytope 190 191 CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 192 */ 193 static inline PetscErrorCode CGNSElementTypeGetTopologyInfo(CGNS_ENUMT(ElementType_t) cellType, DMPolytopeType *dmcelltype, PetscInt *numCorners, PetscInt *dim) 194 { 195 DMPolytopeType _dmcelltype; 196 197 PetscFunctionBeginUser; 198 switch (cellType) { 199 case CGNS_ENUMV(BAR_2): 200 case CGNS_ENUMV(BAR_3): 201 case CGNS_ENUMV(BAR_4): 202 case CGNS_ENUMV(BAR_5): 203 _dmcelltype = DM_POLYTOPE_SEGMENT; 204 break; 205 case CGNS_ENUMV(TRI_3): 206 case CGNS_ENUMV(TRI_6): 207 case CGNS_ENUMV(TRI_9): 208 case CGNS_ENUMV(TRI_10): 209 case CGNS_ENUMV(TRI_12): 210 case CGNS_ENUMV(TRI_15): 211 _dmcelltype = DM_POLYTOPE_TRIANGLE; 212 break; 213 case CGNS_ENUMV(QUAD_4): 214 case CGNS_ENUMV(QUAD_8): 215 case CGNS_ENUMV(QUAD_9): 216 case CGNS_ENUMV(QUAD_12): 217 case CGNS_ENUMV(QUAD_16): 218 case CGNS_ENUMV(QUAD_P4_16): 219 case CGNS_ENUMV(QUAD_25): 220 _dmcelltype = DM_POLYTOPE_QUADRILATERAL; 221 break; 222 case CGNS_ENUMV(TETRA_4): 223 case CGNS_ENUMV(TETRA_10): 224 case CGNS_ENUMV(TETRA_16): 225 case CGNS_ENUMV(TETRA_20): 226 case CGNS_ENUMV(TETRA_22): 227 case CGNS_ENUMV(TETRA_34): 228 case CGNS_ENUMV(TETRA_35): 229 _dmcelltype = DM_POLYTOPE_TETRAHEDRON; 230 break; 231 case CGNS_ENUMV(PYRA_5): 232 case CGNS_ENUMV(PYRA_13): 233 case CGNS_ENUMV(PYRA_14): 234 case CGNS_ENUMV(PYRA_21): 235 case CGNS_ENUMV(PYRA_29): 236 case CGNS_ENUMV(PYRA_P4_29): 237 case CGNS_ENUMV(PYRA_30): 238 case CGNS_ENUMV(PYRA_50): 239 case CGNS_ENUMV(PYRA_55): 240 _dmcelltype = DM_POLYTOPE_PYRAMID; 241 break; 242 case CGNS_ENUMV(PENTA_6): 243 case CGNS_ENUMV(PENTA_15): 244 case CGNS_ENUMV(PENTA_18): 245 case CGNS_ENUMV(PENTA_24): 246 case CGNS_ENUMV(PENTA_33): 247 case CGNS_ENUMV(PENTA_38): 248 case CGNS_ENUMV(PENTA_40): 249 case CGNS_ENUMV(PENTA_66): 250 case CGNS_ENUMV(PENTA_75): 251 _dmcelltype = DM_POLYTOPE_TRI_PRISM; 252 break; 253 case CGNS_ENUMV(HEXA_8): 254 case CGNS_ENUMV(HEXA_20): 255 case CGNS_ENUMV(HEXA_27): 256 case CGNS_ENUMV(HEXA_32): 257 case CGNS_ENUMV(HEXA_44): 258 case CGNS_ENUMV(HEXA_56): 259 case CGNS_ENUMV(HEXA_64): 260 case CGNS_ENUMV(HEXA_98): 261 case CGNS_ENUMV(HEXA_125): 262 _dmcelltype = DM_POLYTOPE_HEXAHEDRON; 263 break; 264 case CGNS_ENUMV(MIXED): 265 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 266 default: 267 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: %d", (int)cellType); 268 } 269 270 if (dmcelltype) *dmcelltype = _dmcelltype; 271 if (numCorners) *numCorners = DMPolytopeTypeGetNumVertices(_dmcelltype); 272 if (dim) *dim = DMPolytopeTypeGetDim(_dmcelltype); 273 PetscFunctionReturn(PETSC_SUCCESS); 274 } 275 276 /* 277 Input Parameters: 278 + cellType - The CGNS-defined cell type 279 280 Output Parameters: 281 + numClosure - Number of nodes that define the function space on the cell 282 - pOrder - The polynomial order of the cell 283 284 CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 285 286 Note: we only support "full" elements, ie. not seredipity elements 287 */ 288 static inline PetscErrorCode CGNSElementTypeGetDiscretizationInfo(CGNS_ENUMT(ElementType_t) cellType, PetscInt *numClosure, PetscInt *pOrder) 289 { 290 PetscInt _numClosure, _pOrder; 291 292 PetscFunctionBeginUser; 293 switch (cellType) { 294 case CGNS_ENUMV(BAR_2): 295 _numClosure = 2; 296 _pOrder = 1; 297 break; 298 case CGNS_ENUMV(BAR_3): 299 _numClosure = 3; 300 _pOrder = 2; 301 break; 302 case CGNS_ENUMV(BAR_4): 303 _numClosure = 4; 304 _pOrder = 3; 305 break; 306 case CGNS_ENUMV(BAR_5): 307 _numClosure = 5; 308 _pOrder = 4; 309 break; 310 case CGNS_ENUMV(TRI_3): 311 _numClosure = 3; 312 _pOrder = 1; 313 break; 314 case CGNS_ENUMV(TRI_6): 315 _numClosure = 6; 316 _pOrder = 2; 317 break; 318 case CGNS_ENUMV(TRI_10): 319 _numClosure = 10; 320 _pOrder = 3; 321 break; 322 case CGNS_ENUMV(TRI_15): 323 _numClosure = 15; 324 _pOrder = 4; 325 break; 326 case CGNS_ENUMV(QUAD_4): 327 _numClosure = 4; 328 _pOrder = 1; 329 break; 330 case CGNS_ENUMV(QUAD_9): 331 _numClosure = 9; 332 _pOrder = 2; 333 break; 334 case CGNS_ENUMV(QUAD_16): 335 _numClosure = 16; 336 _pOrder = 3; 337 break; 338 case CGNS_ENUMV(QUAD_25): 339 _numClosure = 25; 340 _pOrder = 4; 341 break; 342 case CGNS_ENUMV(TETRA_4): 343 _numClosure = 4; 344 _pOrder = 1; 345 break; 346 case CGNS_ENUMV(TETRA_10): 347 _numClosure = 10; 348 _pOrder = 2; 349 break; 350 case CGNS_ENUMV(TETRA_20): 351 _numClosure = 20; 352 _pOrder = 3; 353 break; 354 case CGNS_ENUMV(TETRA_35): 355 _numClosure = 35; 356 _pOrder = 4; 357 break; 358 case CGNS_ENUMV(PYRA_5): 359 _numClosure = 5; 360 _pOrder = 1; 361 break; 362 case CGNS_ENUMV(PYRA_14): 363 _numClosure = 14; 364 _pOrder = 2; 365 break; 366 case CGNS_ENUMV(PYRA_30): 367 _numClosure = 30; 368 _pOrder = 3; 369 break; 370 case CGNS_ENUMV(PYRA_55): 371 _numClosure = 55; 372 _pOrder = 4; 373 break; 374 case CGNS_ENUMV(PENTA_6): 375 _numClosure = 6; 376 _pOrder = 1; 377 break; 378 case CGNS_ENUMV(PENTA_18): 379 _numClosure = 18; 380 _pOrder = 2; 381 break; 382 case CGNS_ENUMV(PENTA_40): 383 _numClosure = 40; 384 _pOrder = 3; 385 break; 386 case CGNS_ENUMV(PENTA_75): 387 _numClosure = 75; 388 _pOrder = 4; 389 break; 390 case CGNS_ENUMV(HEXA_8): 391 _numClosure = 8; 392 _pOrder = 1; 393 break; 394 case CGNS_ENUMV(HEXA_27): 395 _numClosure = 27; 396 _pOrder = 2; 397 break; 398 case CGNS_ENUMV(HEXA_64): 399 _numClosure = 64; 400 _pOrder = 3; 401 break; 402 case CGNS_ENUMV(HEXA_125): 403 _numClosure = 125; 404 _pOrder = 4; 405 break; 406 case CGNS_ENUMV(MIXED): 407 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 408 default: 409 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported or Invalid cell type %d", (int)cellType); 410 } 411 if (numClosure) *numClosure = _numClosure; 412 if (pOrder) *pOrder = _pOrder; 413 PetscFunctionReturn(PETSC_SUCCESS); 414 } 415 416 static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd) 417 { 418 PetscFunctionBegin; 419 switch (pd) { 420 case PETSC_FLOAT: 421 *cd = CGNS_ENUMV(RealSingle); 422 break; 423 case PETSC_DOUBLE: 424 *cd = CGNS_ENUMV(RealDouble); 425 break; 426 case PETSC_COMPLEX: 427 *cd = CGNS_ENUMV(ComplexDouble); 428 break; 429 default: 430 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]); 431 } 432 PetscFunctionReturn(PETSC_SUCCESS); 433 } 434 435 PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 436 { 437 int cgid = -1; 438 PetscBool use_parallel_viewer = PETSC_FALSE; 439 440 PetscFunctionBegin; 441 PetscAssertPointer(filename, 2); 442 PetscCall(PetscViewerCGNSRegisterLogEvents_Internal()); 443 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 444 445 if (use_parallel_viewer) { 446 PetscCallCGNS(cgp_mpi_comm(comm)); 447 PetscCallCGNSOpen(cgp_open(filename, CG_MODE_READ, &cgid), 0, 0); 448 PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cgp_open(\"%s\",...) did not return a valid file ID", filename); 449 PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 450 PetscCallCGNSClose(cgp_close(cgid), 0, 0); 451 } else { 452 PetscCallCGNSOpen(cg_open(filename, CG_MODE_READ, &cgid), 0, 0); 453 PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 454 PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 455 PetscCallCGNSClose(cg_close(cgid), 0, 0); 456 } 457 PetscFunctionReturn(PETSC_SUCCESS); 458 } 459 460 PetscErrorCode DMPlexCreateCGNS_Internal_Serial(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 461 { 462 PetscMPIInt num_proc, rank; 463 DM cdm; 464 DMLabel label; 465 PetscSection coordSection; 466 Vec coordinates; 467 PetscScalar *coords; 468 PetscInt *cellStart, *vertStart, v; 469 PetscInt labelIdRange[2], labelId; 470 /* Read from file */ 471 char basename[CGIO_MAX_NAME_LENGTH + 1]; 472 char buffer[CGIO_MAX_NAME_LENGTH + 1]; 473 int dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0; 474 int nzones = 0; 475 const int B = 1; // Only support single base 476 477 PetscFunctionBegin; 478 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 479 PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 480 PetscCall(DMCreate(comm, dm)); 481 PetscCall(DMSetType(*dm, DMPLEX)); 482 483 /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 484 if (rank == 0) { 485 int nbases, z; 486 487 PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 488 PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 489 PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0); 490 PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0); 491 PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart)); 492 for (z = 1; z <= nzones; ++z) { 493 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 494 495 PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 496 numVertices += sizes[0]; 497 numCells += sizes[1]; 498 cellStart[z] += sizes[1] + cellStart[z - 1]; 499 vertStart[z] += sizes[0] + vertStart[z - 1]; 500 } 501 for (z = 1; z <= nzones; ++z) vertStart[z] += numCells; 502 coordDim = dim; 503 } 504 PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm)); 505 PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 506 PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm)); 507 PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm)); 508 509 PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 510 PetscCall(DMSetDimension(*dm, dim)); 511 PetscCall(DMCreateLabel(*dm, "celltype")); 512 PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 513 514 /* Read zone information */ 515 if (rank == 0) { 516 int z, c, c_loc; 517 518 /* Read the cell set connectivity table and build mesh topology 519 CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 520 /* First set sizes */ 521 for (z = 1, c = 0; z <= nzones; ++z) { 522 CGNS_ENUMT(ZoneType_t) zonetype; 523 int nsections; 524 CGNS_ENUMT(ElementType_t) cellType; 525 cgsize_t start, end; 526 int nbndry, parentFlag; 527 PetscInt numCorners, pOrder; 528 DMPolytopeType ctype; 529 const int S = 1; // Only support single section 530 531 PetscCallCGNSRead(cg_zone_type(cgid, B, z, &zonetype), *dm, 0); 532 PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS"); 533 PetscCallCGNSRead(cg_nsections(cgid, B, z, &nsections), *dm, 0); 534 PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections); 535 PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 536 if (cellType == CGNS_ENUMV(MIXED)) { 537 cgsize_t elementDataSize, *elements; 538 PetscInt off; 539 540 PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 541 PetscCall(PetscMalloc1(elementDataSize, &elements)); 542 PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 543 for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 544 PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[off], &ctype, &numCorners, NULL)); 545 PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[off], NULL, &pOrder)); 546 PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 547 PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 548 PetscCall(DMPlexSetCellType(*dm, c, ctype)); 549 off += numCorners + 1; 550 } 551 PetscCall(PetscFree(elements)); 552 } else { 553 PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &ctype, &numCorners, NULL)); 554 PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 555 PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 556 for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 557 PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 558 PetscCall(DMPlexSetCellType(*dm, c, ctype)); 559 } 560 } 561 } 562 for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 563 } 564 565 PetscCall(DMSetUp(*dm)); 566 567 PetscCall(DMCreateLabel(*dm, "zone")); 568 if (rank == 0) { 569 int z, c, c_loc, v_loc; 570 571 PetscCall(DMGetLabel(*dm, "zone", &label)); 572 for (z = 1, c = 0; z <= nzones; ++z) { 573 CGNS_ENUMT(ElementType_t) cellType; 574 cgsize_t elementDataSize, *elements, start, end; 575 int nbndry, parentFlag; 576 PetscInt *cone, numc, numCorners, maxCorners = 27, pOrder; 577 const int S = 1; // Only support single section 578 579 PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 580 numc = end - start; 581 PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 582 PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone)); 583 PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 584 if (cellType == CGNS_ENUMV(MIXED)) { 585 /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 586 for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 587 PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &numCorners, NULL)); 588 PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &pOrder)); 589 PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 590 ++v; 591 for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 592 PetscCall(DMPlexReorderCell(*dm, c, cone)); 593 PetscCall(DMPlexSetCone(*dm, c, cone)); 594 PetscCall(DMLabelSetValue(label, c, z)); 595 } 596 } else { 597 PetscCall(CGNSElementTypeGetTopologyInfo(cellType, NULL, &numCorners, NULL)); 598 PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 599 PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 600 /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 601 for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 602 for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 603 PetscCall(DMPlexReorderCell(*dm, c, cone)); 604 PetscCall(DMPlexSetCone(*dm, c, cone)); 605 PetscCall(DMLabelSetValue(label, c, z)); 606 } 607 } 608 PetscCall(PetscFree2(elements, cone)); 609 } 610 } 611 612 PetscCall(DMPlexSymmetrize(*dm)); 613 PetscCall(DMPlexStratify(*dm)); 614 if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 615 616 /* Read coordinates */ 617 PetscCall(DMSetCoordinateDim(*dm, coordDim)); 618 PetscCall(DMGetCoordinateDM(*dm, &cdm)); 619 PetscCall(DMGetLocalSection(cdm, &coordSection)); 620 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 621 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim)); 622 PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 623 for (v = numCells; v < numCells + numVertices; ++v) { 624 PetscCall(PetscSectionSetDof(coordSection, v, dim)); 625 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim)); 626 } 627 PetscCall(PetscSectionSetUp(coordSection)); 628 629 PetscCall(DMCreateLocalVector(cdm, &coordinates)); 630 PetscCall(VecGetArray(coordinates, &coords)); 631 if (rank == 0) { 632 PetscInt off = 0; 633 float *x[3]; 634 int z, d; 635 636 PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2])); 637 for (z = 1; z <= nzones; ++z) { 638 CGNS_ENUMT(DataType_t) datatype; 639 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 640 cgsize_t range_min[3] = {1, 1, 1}; 641 cgsize_t range_max[3] = {1, 1, 1}; 642 int ngrids, ncoords; 643 644 PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 645 range_max[0] = sizes[0]; 646 PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0); 647 PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 648 PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0); 649 PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 650 for (d = 0; d < coordDim; ++d) { 651 PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1 + d, &datatype, buffer), *dm, 0); 652 PetscCallCGNSReadData(cg_coord_read(cgid, B, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]), *dm, 0); 653 } 654 if (coordDim >= 1) { 655 for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v]; 656 } 657 if (coordDim >= 2) { 658 for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v]; 659 } 660 if (coordDim >= 3) { 661 for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v]; 662 } 663 off += sizes[0]; 664 } 665 PetscCall(PetscFree3(x[0], x[1], x[2])); 666 } 667 PetscCall(VecRestoreArray(coordinates, &coords)); 668 669 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 670 PetscCall(VecSetBlockSize(coordinates, coordDim)); 671 PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 672 PetscCall(VecDestroy(&coordinates)); 673 674 /* Read boundary conditions */ 675 PetscCall(DMGetNumLabels(*dm, &labelIdRange[0])); 676 if (rank == 0) { 677 CGNS_ENUMT(BCType_t) bctype; 678 CGNS_ENUMT(DataType_t) datatype; 679 CGNS_ENUMT(PointSetType_t) pointtype; 680 cgsize_t *points; 681 PetscReal *normals; 682 int normal[3]; 683 char *bcname = buffer; 684 cgsize_t npoints, nnormals; 685 int z, nbc, bc, c, ndatasets; 686 687 for (z = 1; z <= nzones; ++z) { 688 PetscCallCGNSRead(cg_nbocos(cgid, B, z, &nbc), *dm, 0); 689 for (bc = 1; bc <= nbc; ++bc) { 690 PetscCallCGNSRead(cg_boco_info(cgid, B, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets), *dm, 0); 691 PetscCall(DMCreateLabel(*dm, bcname)); 692 PetscCall(DMGetLabel(*dm, bcname, &label)); 693 PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals)); 694 PetscCallCGNSReadData(cg_boco_read(cgid, B, z, bc, points, (void *)normals), *dm, 0); 695 if (pointtype == CGNS_ENUMV(ElementRange)) { 696 // Range of cells: assuming half-open interval 697 for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 698 } else if (pointtype == CGNS_ENUMV(ElementList)) { 699 // List of cells 700 for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 701 } else if (pointtype == CGNS_ENUMV(PointRange)) { 702 CGNS_ENUMT(GridLocation_t) gridloc; 703 704 // List of points: 705 PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 706 PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 707 // Range of points: assuming half-open interval 708 for (c = points[0]; c < points[1]; ++c) { 709 if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1)); 710 else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 711 } 712 } else if (pointtype == CGNS_ENUMV(PointList)) { 713 CGNS_ENUMT(GridLocation_t) gridloc; 714 715 // List of points: 716 PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 717 PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 718 for (c = 0; c < npoints; ++c) { 719 if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1)); 720 else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 721 } 722 } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype); 723 PetscCall(PetscFree2(points, normals)); 724 } 725 } 726 PetscCall(PetscFree2(cellStart, vertStart)); 727 } 728 PetscCall(DMGetNumLabels(*dm, &labelIdRange[1])); 729 PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm)); 730 731 /* Create BC labels at all processes */ 732 for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 733 char *labelName = buffer; 734 size_t len = sizeof(buffer); 735 const char *locName; 736 737 if (rank == 0) { 738 PetscCall(DMGetLabelByNum(*dm, labelId, &label)); 739 PetscCall(PetscObjectGetName((PetscObject)label, &locName)); 740 PetscCall(PetscStrncpy(labelName, locName, len)); 741 } 742 PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm)); 743 PetscCallMPI(DMCreateLabel(*dm, labelName)); 744 } 745 PetscFunctionReturn(PETSC_SUCCESS); 746 } 747 748 PetscErrorCode DMPlexCreateCGNS_Internal_Parallel(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 749 { 750 PetscMPIInt num_proc, rank; 751 /* Read from file */ 752 char basename[CGIO_MAX_NAME_LENGTH + 1]; 753 char buffer[CGIO_MAX_NAME_LENGTH + 1]; 754 int dim = 0, physDim = 0, coordDim = 0; 755 PetscInt NVertices = 0, NCells = 0; 756 int nzones = 0, nbases; 757 int z = 1; // Only supports single zone files 758 int B = 1; // Only supports single base 759 760 PetscFunctionBegin; 761 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 762 PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 763 PetscCall(DMCreate(comm, dm)); 764 PetscCall(DMSetType(*dm, DMPLEX)); 765 766 PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 767 PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 768 // From the CGNS web page cell_dim phys_dim (embedding space in PETSc) CGNS defines as length of spatial vectors/components) 769 PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0); 770 PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0); 771 PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Parallel reader limited to one zone, not %d", nzones); 772 { 773 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 774 775 PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 776 NVertices = sizes[0]; 777 NCells = sizes[1]; 778 } 779 780 PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 781 PetscCall(DMSetDimension(*dm, dim)); 782 coordDim = dim; 783 784 // This is going to be a headache for mixed-topology and multiple sections. We may have to restore reading the data twice (once before the SetChart 785 // call to get this right but continuing for now with single section, single topology, one zone. 786 // establish element ranges for my rank 787 PetscInt mystarte, myende, mystartv, myendv, myownede, myownedv; 788 PetscLayout elem_map, vtx_map; 789 PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NCells, 1, &elem_map)); 790 PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 791 PetscCall(PetscLayoutGetRange(elem_map, &mystarte, &myende)); 792 PetscCall(PetscLayoutGetLocalSize(elem_map, &myownede)); 793 PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 794 PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 795 796 // -- Build Plex in parallel 797 DMPolytopeType dm_cell_type = DM_POLYTOPE_UNKNOWN; 798 PetscInt pOrder = 1, numClosure = -1; 799 cgsize_t *elements; 800 { 801 int nsections; 802 PetscInt *elementsQ1, numCorners = -1; 803 const int *perm; 804 cgsize_t start, end; // Throwaway 805 806 cg_nsections(cgid, B, z, &nsections); 807 // Read element connectivity 808 for (int index_sect = 1; index_sect <= nsections; index_sect++) { 809 int nbndry, parentFlag; 810 PetscInt cell_dim; 811 CGNS_ENUMT(ElementType_t) cellType; 812 813 PetscCallCGNSRead(cg_section_read(cgid, B, z, index_sect, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 814 815 PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &dm_cell_type, &numCorners, &cell_dim)); 816 // Skip over element that are not max dimension (ie. boundary elements) 817 if (cell_dim != dim) continue; 818 PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, &numClosure, &pOrder)); 819 PetscCall(PetscMalloc1(myownede * numClosure, &elements)); 820 PetscCallCGNSReadData(cgp_elements_read_data(cgid, B, z, index_sect, mystarte + 1, myende, elements), *dm, 0); 821 for (PetscInt v = 0; v < myownede * numClosure; ++v) elements[v] -= 1; // 0 based 822 break; 823 } 824 825 // Create corners-only connectivity 826 PetscCall(PetscMalloc1(myownede * numCorners, &elementsQ1)); 827 PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numCorners, NULL, &perm)); 828 for (PetscInt e = 0; e < myownede; ++e) { 829 for (PetscInt v = 0; v < numCorners; ++v) elementsQ1[e * numCorners + perm[v]] = elements[e * numClosure + v]; 830 } 831 832 // Build cell-vertex Plex 833 PetscCall(DMPlexBuildFromCellListParallel(*dm, myownede, myownedv, NVertices, numCorners, elementsQ1, NULL, NULL)); 834 PetscCall(DMViewFromOptions(*dm, NULL, "-corner_dm_view")); 835 PetscCall(PetscFree(elementsQ1)); 836 } 837 838 if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 839 840 // -- Create SF for naive nodal-data read to elements 841 PetscSF plex_to_cgns_sf; 842 { 843 PetscInt nleaves, num_comp; 844 PetscInt *leaf, num_leaves = 0; 845 PetscInt cStart, cEnd; 846 const int *perm; 847 PetscSF cgns_to_local_sf; 848 PetscSection local_section; 849 PetscFE fe; 850 851 // sfNatural requires PetscSection to handle DMDistribute, so we use PetscFE to define the section 852 // Use number of components = 1 to work with just the nodes themselves 853 PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, dm_cell_type, pOrder, PETSC_DETERMINE, &fe)); 854 PetscCall(PetscObjectSetName((PetscObject)fe, "FE for sfNatural")); 855 PetscCall(DMAddField(*dm, NULL, (PetscObject)fe)); 856 PetscCall(DMCreateDS(*dm)); 857 PetscCall(PetscFEDestroy(&fe)); 858 859 PetscCall(DMGetLocalSection(*dm, &local_section)); 860 PetscCall(PetscSectionViewFromOptions(local_section, NULL, "-fe_natural_section_view")); 861 PetscCall(PetscSectionGetFieldComponents(local_section, 0, &num_comp)); 862 PetscCall(PetscSectionGetStorageSize(local_section, &nleaves)); 863 nleaves /= num_comp; 864 PetscCall(PetscMalloc1(nleaves, &leaf)); 865 for (PetscInt i = 0; i < nleaves; i++) leaf[i] = -1; 866 867 // Get the permutation from CGNS closure numbering to PLEX closure numbering 868 PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numClosure, NULL, &perm)); 869 PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 870 for (PetscInt cell = cStart; cell < cEnd; ++cell) { 871 PetscInt num_closure_dof, *closure_idx = NULL; 872 873 PetscCall(DMPlexGetClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 874 PetscAssert(numClosure * num_comp == num_closure_dof, comm, PETSC_ERR_PLIB, "Closure dof size does not match polytope"); 875 for (PetscInt i = 0; i < numClosure; i++) { 876 PetscInt li = closure_idx[perm[i] * num_comp] / num_comp; 877 if (li < 0) continue; 878 879 PetscInt cgns_idx = elements[cell * numClosure + i]; 880 if (leaf[li] == -1) { 881 leaf[li] = cgns_idx; 882 num_leaves++; 883 } else PetscAssert(leaf[li] == cgns_idx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf does not match previously set"); 884 } 885 PetscCall(DMPlexRestoreClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 886 } 887 PetscAssert(num_leaves == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf count in closure does not match nleaves"); 888 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)*dm), &cgns_to_local_sf)); 889 PetscCall(PetscSFSetGraphLayout(cgns_to_local_sf, vtx_map, nleaves, NULL, PETSC_USE_POINTER, leaf)); 890 PetscCall(PetscObjectSetName((PetscObject)cgns_to_local_sf, "CGNS to Plex SF")); 891 PetscCall(PetscSFViewFromOptions(cgns_to_local_sf, NULL, "-CGNStoPlex_sf_view")); 892 PetscCall(PetscFree(leaf)); 893 PetscCall(PetscFree(elements)); 894 895 { // Convert cgns_to_local to global_to_cgns 896 PetscSF sectionsf, cgns_to_global_sf; 897 898 PetscCall(DMGetSectionSF(*dm, §ionsf)); 899 PetscCall(PetscSFComposeInverse(cgns_to_local_sf, sectionsf, &cgns_to_global_sf)); 900 PetscCall(PetscSFDestroy(&cgns_to_local_sf)); 901 PetscCall(PetscSFCreateInverseSF(cgns_to_global_sf, &plex_to_cgns_sf)); 902 PetscCall(PetscObjectSetName((PetscObject)plex_to_cgns_sf, "Global Plex to CGNS")); 903 PetscCall(PetscSFDestroy(&cgns_to_global_sf)); 904 } 905 } 906 907 { // -- Set coordinates for DM 908 PetscScalar *coords; 909 float *x[3]; 910 double *xd[3]; 911 PetscBool read_with_double; 912 PetscFE cfe; 913 914 // Setup coordinate space first. Use pOrder here for isoparametric; revisit with CPEX-0045 High Order. 915 PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, coordDim, dm_cell_type, pOrder, PETSC_DETERMINE, &cfe)); 916 PetscCall(DMSetCoordinateDisc(*dm, cfe, PETSC_FALSE, PETSC_FALSE)); 917 PetscCall(PetscFEDestroy(&cfe)); 918 919 { // Determine if coords are written in single or double precision 920 CGNS_ENUMT(DataType_t) datatype; 921 922 PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1, &datatype, buffer), *dm, 0); 923 read_with_double = datatype == CGNS_ENUMV(RealDouble) ? PETSC_TRUE : PETSC_FALSE; 924 } 925 926 // Read coords from file and set into component-major ordering 927 if (read_with_double) PetscCall(PetscMalloc3(myownedv, &xd[0], myownedv, &xd[1], myownedv, &xd[2])); 928 else PetscCall(PetscMalloc3(myownedv, &x[0], myownedv, &x[1], myownedv, &x[2])); 929 PetscCall(PetscMalloc1(myownedv * coordDim, &coords)); 930 { 931 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 932 cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 933 cgsize_t range_max[3] = {myendv, 1, 1}; 934 int ngrids, ncoords; 935 936 PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 937 PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0); 938 PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 939 PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0); 940 PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 941 if (read_with_double) { 942 for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, B, z, (d + 1), range_min, range_max, xd[d]), *dm, 0); 943 if (coordDim >= 1) { 944 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = xd[0][v]; 945 } 946 if (coordDim >= 2) { 947 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = xd[1][v]; 948 } 949 if (coordDim >= 3) { 950 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = xd[2][v]; 951 } 952 } else { 953 for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, 1, z, (d + 1), range_min, range_max, x[d]), *dm, 0); 954 if (coordDim >= 1) { 955 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = x[0][v]; 956 } 957 if (coordDim >= 2) { 958 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = x[1][v]; 959 } 960 if (coordDim >= 3) { 961 for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = x[2][v]; 962 } 963 } 964 } 965 if (read_with_double) PetscCall(PetscFree3(xd[0], xd[1], xd[2])); 966 else PetscCall(PetscFree3(x[0], x[1], x[2])); 967 968 { // Reduce CGNS-ordered coordinate nodes to Plex ordering, and set DM's coordinates 969 Vec coord_global; 970 MPI_Datatype unit; 971 PetscScalar *coord_global_array; 972 DM cdm; 973 974 PetscCall(DMGetCoordinateDM(*dm, &cdm)); 975 PetscCall(DMCreateGlobalVector(cdm, &coord_global)); 976 PetscCall(VecGetArrayWrite(coord_global, &coord_global_array)); 977 PetscCallMPI(MPI_Type_contiguous(coordDim, MPIU_SCALAR, &unit)); 978 PetscCallMPI(MPI_Type_commit(&unit)); 979 PetscCall(PetscSFReduceBegin(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 980 PetscCall(PetscSFReduceEnd(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 981 PetscCall(VecRestoreArrayWrite(coord_global, &coord_global_array)); 982 PetscCallMPI(MPI_Type_free(&unit)); 983 PetscCall(DMSetCoordinates(*dm, coord_global)); 984 PetscCall(VecDestroy(&coord_global)); 985 } 986 PetscCall(PetscFree(coords)); 987 } 988 989 // -- Set sfNatural for solution vectors in CGNS file 990 // NOTE: We set sfNatural to be the map between the original CGNS ordering of nodes and the Plex ordering of nodes. 991 PetscCall(PetscSFViewFromOptions(plex_to_cgns_sf, NULL, "-sfNatural_init_view")); 992 PetscCall(DMSetNaturalSF(*dm, plex_to_cgns_sf)); 993 PetscCall(DMSetUseNatural(*dm, PETSC_TRUE)); 994 PetscCall(PetscSFDestroy(&plex_to_cgns_sf)); 995 996 PetscCall(PetscLayoutDestroy(&elem_map)); 997 PetscCall(PetscLayoutDestroy(&vtx_map)); 998 PetscFunctionReturn(PETSC_SUCCESS); 999 } 1000 1001 // node_l2g must be freed 1002 static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g) 1003 { 1004 PetscSection local_section; 1005 PetscSF point_sf; 1006 PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes; 1007 PetscMPIInt comm_size; 1008 const PetscInt *ilocal, field = 0; 1009 1010 PetscFunctionBegin; 1011 *num_local_nodes = -1; 1012 *num_global_nodes = -1; 1013 *nStart = -1; 1014 *nEnd = -1; 1015 *node_l2g = NULL; 1016 1017 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size)); 1018 PetscCall(DMGetLocalSection(dm, &local_section)); 1019 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1020 PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd)); 1021 PetscCall(DMGetPointSF(dm, &point_sf)); 1022 if (comm_size == 1) nleaves = 0; 1023 else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL)); 1024 PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp)); 1025 1026 PetscInt local_node = 0, owned_node = 0, owned_start = 0; 1027 for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 1028 PetscInt dof; 1029 PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 1030 PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp); 1031 local_node += dof / ncomp; 1032 if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 1033 leaf++; 1034 } else { 1035 owned_node += dof / ncomp; 1036 } 1037 } 1038 PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 1039 PetscCall(PetscMalloc1(pEnd - pStart, &points)); 1040 owned_node = 0; 1041 for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 1042 if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 1043 points[p - pStart] = -1; 1044 leaf++; 1045 continue; 1046 } 1047 PetscInt dof, offset; 1048 PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 1049 PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 1050 PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp); 1051 points[p - pStart] = owned_start + owned_node; 1052 owned_node += dof / ncomp; 1053 } 1054 if (comm_size > 1) { 1055 PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 1056 PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 1057 } 1058 1059 // Set up global indices for each local node 1060 PetscCall(PetscMalloc1(local_node, &nodes)); 1061 for (PetscInt p = spStart; p < spEnd; p++) { 1062 PetscInt dof, offset; 1063 PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 1064 PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 1065 for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n; 1066 } 1067 PetscCall(PetscFree(points)); 1068 *num_local_nodes = local_node; 1069 *nStart = owned_start; 1070 *nEnd = owned_start + owned_node; 1071 PetscCallMPI(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 1072 *node_l2g = nodes; 1073 PetscFunctionReturn(PETSC_SUCCESS); 1074 } 1075 1076 PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer) 1077 { 1078 PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 1079 PetscInt fvGhostStart; 1080 PetscInt topo_dim, coord_dim, num_global_elems; 1081 PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd; 1082 const PetscInt *node_l2g; 1083 Vec coord; 1084 DM colloc_dm, cdm; 1085 PetscMPIInt size; 1086 const char *dm_name; 1087 int base, zone; 1088 cgsize_t isize[3]; 1089 1090 PetscFunctionBegin; 1091 if (!cgv->file_num) { 1092 PetscInt time_step; 1093 PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL)); 1094 PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step)); 1095 } 1096 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 1097 PetscCall(DMGetDimension(dm, &topo_dim)); 1098 PetscCall(DMGetCoordinateDim(dm, &coord_dim)); 1099 PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name)); 1100 PetscCallCGNSWrite(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base), dm, viewer); 1101 PetscCallCGNS(cg_goto(cgv->file_num, base, NULL)); 1102 PetscCallCGNSWrite(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)), dm, viewer); 1103 1104 { 1105 PetscFE fe, fe_coord; 1106 PetscClassId ds_id; 1107 PetscDualSpace dual_space, dual_space_coord; 1108 PetscInt num_fields, field_order = -1, field_order_coord; 1109 PetscBool is_simplex; 1110 PetscCall(DMGetNumFields(dm, &num_fields)); 1111 if (num_fields > 0) { 1112 PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 1113 PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id)); 1114 if (ds_id != PETSCFE_CLASSID) { 1115 fe = NULL; 1116 if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter 1117 else field_order = 1; // assume vertex-based linear elements 1118 } 1119 } else fe = NULL; 1120 if (fe) { 1121 PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 1122 PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order)); 1123 } 1124 PetscCall(DMGetCoordinateDM(dm, &cdm)); 1125 PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord)); 1126 { 1127 PetscClassId id; 1128 PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id)); 1129 if (id != PETSCFE_CLASSID) fe_coord = NULL; 1130 } 1131 if (fe_coord) { 1132 PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord)); 1133 PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord)); 1134 } else field_order_coord = 1; 1135 if (field_order > 0 && field_order != field_order_coord) { 1136 PetscInt quadrature_order = field_order; 1137 PetscCall(DMClone(dm, &colloc_dm)); 1138 { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied 1139 const PetscSF *face_sfs; 1140 PetscInt num_face_sfs; 1141 PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs)); 1142 PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs)); 1143 if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal; 1144 } 1145 PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 1146 PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe)); 1147 PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_FALSE, PETSC_TRUE)); 1148 PetscCall(PetscFEDestroy(&fe)); 1149 } else { 1150 PetscCall(PetscObjectReference((PetscObject)dm)); 1151 colloc_dm = dm; 1152 } 1153 } 1154 PetscCall(DMGetCoordinateDM(colloc_dm, &cdm)); 1155 PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g)); 1156 PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord)); 1157 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1158 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 1159 if (fvGhostStart >= 0) cEnd = fvGhostStart; 1160 num_global_elems = cEnd - cStart; 1161 PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 1162 isize[0] = num_global_nodes; 1163 isize[1] = num_global_elems; 1164 isize[2] = 0; 1165 PetscCallCGNSWrite(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone), dm, viewer); 1166 1167 { 1168 const PetscScalar *X; 1169 PetscScalar *x; 1170 int coord_ids[3]; 1171 1172 PetscCall(VecGetArrayRead(coord, &X)); 1173 for (PetscInt d = 0; d < coord_dim; d++) { 1174 const double exponents[] = {0, 1, 0, 0, 0}; 1175 char coord_name[64]; 1176 PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d)); 1177 PetscCallCGNSWrite(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]), dm, viewer); 1178 PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL)); 1179 PetscCallCGNSWrite(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents), dm, viewer); 1180 } 1181 1182 int section; 1183 cgsize_t e_owned, e_global, e_start, *conn = NULL; 1184 const int *perm; 1185 CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull); 1186 { 1187 PetscCall(PetscMalloc1(nEnd - nStart, &x)); 1188 for (PetscInt d = 0; d < coord_dim; d++) { 1189 for (PetscInt n = 0; n < num_local_nodes; n++) { 1190 PetscInt gn = node_l2g[n]; 1191 if (gn < nStart || nEnd <= gn) continue; 1192 x[gn - nStart] = X[n * coord_dim + d]; 1193 } 1194 // CGNS nodes use 1-based indexing 1195 cgsize_t start = nStart + 1, end = nEnd; 1196 PetscCallCGNSWriteData(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x), dm, viewer); 1197 } 1198 PetscCall(PetscFree(x)); 1199 PetscCall(VecRestoreArrayRead(coord, &X)); 1200 } 1201 1202 e_owned = cEnd - cStart; 1203 if (e_owned > 0) { 1204 DMPolytopeType cell_type; 1205 1206 PetscCall(DMPlexGetCellType(dm, cStart, &cell_type)); 1207 for (PetscInt i = cStart, c = 0; i < cEnd; i++) { 1208 PetscInt closure_dof, *closure_indices, elem_size; 1209 1210 PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 1211 elem_size = closure_dof / coord_dim; 1212 if (!conn) PetscCall(PetscMalloc1(e_owned * elem_size, &conn)); 1213 PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm)); 1214 for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1; 1215 PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 1216 } 1217 } 1218 1219 { // Get global element_type (for ranks that do not have owned elements) 1220 PetscInt local_element_type, global_element_type; 1221 1222 local_element_type = e_owned > 0 ? element_type : -1; 1223 PetscCallMPI(MPIU_Allreduce(&local_element_type, &global_element_type, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 1224 if (local_element_type != -1) PetscCheck(local_element_type == global_element_type, PETSC_COMM_SELF, PETSC_ERR_SUP, "Ranks with different element types not supported"); 1225 element_type = (CGNS_ENUMT(ElementType_t))global_element_type; 1226 } 1227 PetscCallMPI(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 1228 PetscCheck(e_global == num_global_elems, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unexpected number of elements %" PRIdCGSIZE " vs %" PetscInt_FMT, e_global, num_global_elems); 1229 e_start = 0; 1230 PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 1231 PetscCallCGNSWrite(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, §ion), dm, viewer); 1232 PetscCallCGNSWriteData(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn), dm, viewer); 1233 PetscCall(PetscFree(conn)); 1234 1235 cgv->base = base; 1236 cgv->zone = zone; 1237 cgv->node_l2g = node_l2g; 1238 cgv->num_local_nodes = num_local_nodes; 1239 cgv->nStart = nStart; 1240 cgv->nEnd = nEnd; 1241 cgv->eStart = e_start; 1242 cgv->eEnd = e_start + e_owned; 1243 if (1) { 1244 PetscMPIInt rank; 1245 int *efield; 1246 int sol, field; 1247 DMLabel label; 1248 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 1249 PetscCall(PetscMalloc1(e_owned, &efield)); 1250 for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank; 1251 PetscCallCGNSWrite(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol), dm, viewer); 1252 PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field), dm, viewer); 1253 cgsize_t start = e_start + 1, end = e_start + e_owned; 1254 PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 1255 PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 1256 if (label) { 1257 for (PetscInt c = cStart; c < cEnd; c++) { 1258 PetscInt value; 1259 PetscCall(DMLabelGetValue(label, c, &value)); 1260 efield[c - cStart] = value; 1261 } 1262 PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field), dm, viewer); 1263 PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 1264 } 1265 PetscCall(PetscFree(efield)); 1266 } 1267 } 1268 PetscCall(DMDestroy(&colloc_dm)); 1269 PetscFunctionReturn(PETSC_SUCCESS); 1270 } 1271 1272 PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer) 1273 { 1274 PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 1275 DM dm; 1276 PetscSection section; 1277 PetscInt time_step, num_fields, pStart, pEnd, fvGhostStart; 1278 PetscReal time, *time_slot; 1279 size_t *step_slot; 1280 const PetscScalar *v; 1281 char solution_name[PETSC_MAX_PATH_LEN]; 1282 int sol; 1283 1284 PetscFunctionBegin; 1285 PetscCall(VecGetDM(V, &dm)); 1286 PetscCall(DMGetLocalSection(dm, §ion)); 1287 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 1288 PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 1289 if (fvGhostStart >= 0) pEnd = fvGhostStart; 1290 1291 if (!cgv->node_l2g) PetscCall(DMView(dm, viewer)); 1292 if (!cgv->grid_loc) { // Determine if writing to cell-centers or to nodes 1293 PetscInt cStart, cEnd; 1294 PetscInt local_grid_loc, global_grid_loc; 1295 1296 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 1297 if (fvGhostStart >= 0) cEnd = fvGhostStart; 1298 if (cgv->num_local_nodes == 0) local_grid_loc = -1; 1299 else if (cStart == pStart && cEnd == pEnd) local_grid_loc = CGNS_ENUMV(CellCenter); 1300 else local_grid_loc = CGNS_ENUMV(Vertex); 1301 1302 PetscCallMPI(MPIU_Allreduce(&local_grid_loc, &global_grid_loc, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 1303 if (local_grid_loc != -1) 1304 PetscCheck(local_grid_loc == global_grid_loc, PETSC_COMM_SELF, PETSC_ERR_SUP, "Ranks with different grid locations not supported. Local has %" PetscInt_FMT ", allreduce returned %" PetscInt_FMT, local_grid_loc, global_grid_loc); 1305 PetscCheck((global_grid_loc == CGNS_ENUMV(CellCenter)) || (global_grid_loc == CGNS_ENUMV(Vertex)), PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Grid location should only be CellCenter (%d) or Vertex(%d), but have %" PetscInt_FMT, CGNS_ENUMV(CellCenter), CGNS_ENUMV(Vertex), global_grid_loc); 1306 cgv->grid_loc = (CGNS_ENUMT(GridLocation_t))global_grid_loc; 1307 } 1308 if (!cgv->nodal_field) { 1309 switch (cgv->grid_loc) { 1310 case CGNS_ENUMV(Vertex): { 1311 PetscCall(PetscMalloc1(cgv->nEnd - cgv->nStart, &cgv->nodal_field)); 1312 } break; 1313 case CGNS_ENUMV(CellCenter): { 1314 PetscCall(PetscMalloc1(cgv->eEnd - cgv->eStart, &cgv->nodal_field)); 1315 } break; 1316 default: 1317 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 1318 } 1319 } 1320 if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times)); 1321 if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps)); 1322 1323 PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time)); 1324 if (time_step < 0) { 1325 time_step = 0; 1326 time = 0.; 1327 } 1328 // Avoid "Duplicate child name found" error by not writing an already-written solution. 1329 // This usually occurs when a solution is written and then diverges on the very next timestep. 1330 if (time_step == cgv->previous_output_step) PetscFunctionReturn(PETSC_SUCCESS); 1331 1332 PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot)); 1333 *time_slot = time; 1334 PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot)); 1335 *step_slot = cgv->previous_output_step = time_step; 1336 PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step)); 1337 PetscCallCGNSWrite(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, cgv->grid_loc, &sol), V, viewer); 1338 PetscCall(VecGetArrayRead(V, &v)); 1339 PetscCall(PetscSectionGetNumFields(section, &num_fields)); 1340 for (PetscInt field = 0; field < num_fields; field++) { 1341 PetscInt ncomp; 1342 const char *field_name; 1343 PetscCall(PetscSectionGetFieldName(section, field, &field_name)); 1344 PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp)); 1345 for (PetscInt comp = 0; comp < ncomp; comp++) { 1346 int cgfield; 1347 const char *comp_name; 1348 char cgns_field_name[32]; // CGNS max field name is 32 1349 CGNS_ENUMT(DataType_t) datatype; 1350 PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name)); 1351 if (ncomp == 1 && comp_name[0] == '0' && comp_name[1] == '\0' && field_name[0] != '\0') PetscCall(PetscStrncpy(cgns_field_name, field_name, sizeof cgns_field_name)); 1352 else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name)); 1353 else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name)); 1354 PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype)); 1355 PetscCallCGNSWrite(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield), V, viewer); 1356 for (PetscInt p = pStart, n = 0; p < pEnd; p++) { 1357 PetscInt off, dof; 1358 PetscCall(PetscSectionGetFieldDof(section, p, field, &dof)); 1359 if (dof == 0) continue; 1360 PetscCall(PetscSectionGetFieldOffset(section, p, field, &off)); 1361 for (PetscInt c = comp; c < dof; c += ncomp, n++) { 1362 switch (cgv->grid_loc) { 1363 case CGNS_ENUMV(Vertex): { 1364 PetscInt gn = cgv->node_l2g[n]; 1365 if (gn < cgv->nStart || cgv->nEnd <= gn) continue; 1366 cgv->nodal_field[gn - cgv->nStart] = v[off + c]; 1367 } break; 1368 case CGNS_ENUMV(CellCenter): { 1369 cgv->nodal_field[n] = v[off + c]; 1370 } break; 1371 default: 1372 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations"); 1373 } 1374 } 1375 } 1376 // CGNS nodes use 1-based indexing 1377 cgsize_t start, end; 1378 switch (cgv->grid_loc) { 1379 case CGNS_ENUMV(Vertex): { 1380 start = cgv->nStart + 1; 1381 end = cgv->nEnd; 1382 } break; 1383 case CGNS_ENUMV(CellCenter): { 1384 start = cgv->eStart + 1; 1385 end = cgv->eEnd; 1386 } break; 1387 default: 1388 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 1389 } 1390 PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field), V, viewer); 1391 } 1392 } 1393 PetscCall(VecRestoreArrayRead(V, &v)); 1394 PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer)); 1395 PetscFunctionReturn(PETSC_SUCCESS); 1396 } 1397 1398 PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec V, PetscViewer viewer) 1399 { 1400 MPI_Comm comm; 1401 char buffer[CGIO_MAX_NAME_LENGTH + 1]; 1402 PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 1403 int cgid = cgv->file_num; 1404 PetscBool use_parallel_viewer = PETSC_FALSE; 1405 int z = 1; // Only support one zone 1406 int B = 1; // Only support one base 1407 int numComp; 1408 PetscInt V_numComps, mystartv, myendv, myownedv; 1409 1410 PetscFunctionBeginUser; 1411 PetscCall(PetscObjectGetComm((PetscObject)V, &comm)); 1412 1413 PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 1414 PetscCheck(use_parallel_viewer, comm, PETSC_ERR_USER_INPUT, "Cannot use VecLoad with CGNS file in serial reader; use -dm_plex_cgns_parallel to enable parallel reader"); 1415 1416 { // Get CGNS node ownership information 1417 int nbases, nzones; 1418 PetscInt NVertices; 1419 PetscLayout vtx_map; 1420 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 1421 1422 PetscCallCGNSRead(cg_nbases(cgid, &nbases), V, viewer); 1423 PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 1424 PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), V, viewer); 1425 PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "limited to one zone %d", (int)nzones); 1426 1427 PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), V, viewer); 1428 NVertices = sizes[0]; 1429 1430 PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 1431 PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 1432 PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 1433 PetscCall(PetscLayoutDestroy(&vtx_map)); 1434 } 1435 1436 { // -- Read data from file into Vec 1437 PetscScalar *fields = NULL; 1438 PetscSF sfNatural; 1439 1440 { // Check compatibility between sfNatural and the data source and sink 1441 DM dm; 1442 PetscInt nleaves, nroots, V_local_size; 1443 1444 PetscCall(VecGetDM(V, &dm)); 1445 PetscCall(DMGetNaturalSF(dm, &sfNatural)); 1446 PetscCheck(sfNatural, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM of Vec must have sfNatural"); 1447 PetscCall(PetscSFGetGraph(sfNatural, &nroots, &nleaves, NULL, NULL)); 1448 PetscCall(VecGetLocalSize(V, &V_local_size)); 1449 PetscCheck(nleaves == myownedv, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Number of locally owned vertices (% " PetscInt_FMT ") must match number of leaves in sfNatural (% " PetscInt_FMT ")", myownedv, nleaves); 1450 PetscCheck(V_local_size % nroots == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Local Vec size (% " PetscInt_FMT ") not evenly divisible by number of roots in sfNatural (% " PetscInt_FMT ")", V_local_size, nroots); 1451 V_numComps = V_local_size / nroots; 1452 } 1453 1454 { // Read data into component-major ordering 1455 int isol, numSols; 1456 CGNS_ENUMT(DataType_t) datatype; 1457 double *fields_CGNS; 1458 1459 PetscCallCGNSRead(cg_nsols(cgid, B, z, &numSols), V, viewer); 1460 PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &isol)); 1461 PetscCallCGNSRead(cg_nfields(cgid, B, z, isol, &numComp), V, viewer); 1462 PetscCheck(V_numComps == numComp, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Vec sized for % " PetscInt_FMT " components per node, but file has %d components per node", V_numComps, numComp); 1463 1464 cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 1465 cgsize_t range_max[3] = {myendv, 1, 1}; 1466 PetscCall(PetscMalloc1(myownedv * numComp, &fields_CGNS)); 1467 PetscCall(PetscMalloc1(myownedv * numComp, &fields)); 1468 for (int d = 0; d < numComp; ++d) { 1469 PetscCallCGNSRead(cg_field_info(cgid, B, z, isol, (d + 1), &datatype, buffer), V, viewer); 1470 PetscCheck(datatype == CGNS_ENUMV(RealDouble), PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Field %s in file is not of type double", buffer); 1471 PetscCallCGNSReadData(cgp_field_read_data(cgid, B, z, isol, (d + 1), range_min, range_max, &fields_CGNS[d * myownedv]), V, viewer); 1472 } 1473 for (int d = 0; d < numComp; ++d) { 1474 for (PetscInt v = 0; v < myownedv; ++v) fields[v * numComp + d] = fields_CGNS[d * myownedv + v]; 1475 } 1476 PetscCall(PetscFree(fields_CGNS)); 1477 } 1478 1479 { // Reduce fields into Vec array 1480 PetscScalar *V_array; 1481 MPI_Datatype fieldtype; 1482 1483 PetscCall(VecGetArrayWrite(V, &V_array)); 1484 PetscCallMPI(MPI_Type_contiguous(numComp, MPIU_SCALAR, &fieldtype)); 1485 PetscCallMPI(MPI_Type_commit(&fieldtype)); 1486 PetscCall(PetscSFReduceBegin(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1487 PetscCall(PetscSFReduceEnd(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1488 PetscCallMPI(MPI_Type_free(&fieldtype)); 1489 PetscCall(VecRestoreArrayWrite(V, &V_array)); 1490 } 1491 PetscCall(PetscFree(fields)); 1492 } 1493 PetscFunctionReturn(PETSC_SUCCESS); 1494 } 1495