15f34f2dcSJed Brown #define PETSCDM_DLL 25f34f2dcSJed Brown #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 35f34f2dcSJed Brown #include <petsc/private/viewercgnsimpl.h> 45f34f2dcSJed Brown 55f34f2dcSJed Brown #include <pcgnslib.h> 65f34f2dcSJed Brown #include <cgns_io.h> 75f34f2dcSJed Brown 85f34f2dcSJed Brown #if !defined(CGNS_ENUMT) 95f34f2dcSJed Brown #define CGNS_ENUMT(a) a 105f34f2dcSJed Brown #endif 115f34f2dcSJed Brown #if !defined(CGNS_ENUMV) 125f34f2dcSJed Brown #define CGNS_ENUMV(a) a 135f34f2dcSJed Brown #endif 145f34f2dcSJed Brown // Permute plex closure ordering to CGNS 15d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm) 16d71ae5a4SJacob Faibussowitsch { 17472b9844SJames Wright CGNS_ENUMT(ElementType_t) element_type_tmp; 18472b9844SJames Wright 19472b9844SJames Wright // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 205f34f2dcSJed Brown static const int bar_2[2] = {0, 1}; 215f34f2dcSJed Brown static const int bar_3[3] = {1, 2, 0}; 225f34f2dcSJed Brown static const int bar_4[4] = {2, 3, 0, 1}; 235f34f2dcSJed Brown static const int bar_5[5] = {3, 4, 0, 1, 2}; 245f34f2dcSJed Brown static const int tri_3[3] = {0, 1, 2}; 255f34f2dcSJed Brown static const int tri_6[6] = {3, 4, 5, 0, 1, 2}; 265f34f2dcSJed Brown static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0}; 275f34f2dcSJed Brown static const int quad_4[4] = {0, 1, 2, 3}; 285f34f2dcSJed Brown static const int quad_9[9] = { 295f34f2dcSJed Brown 5, 6, 7, 8, // vertices 305f34f2dcSJed Brown 1, 2, 3, 4, // edges 315f34f2dcSJed Brown 0, // center 325f34f2dcSJed Brown }; 335f34f2dcSJed Brown static const int quad_16[] = { 345f34f2dcSJed Brown 12, 13, 14, 15, // vertices 355f34f2dcSJed Brown 4, 5, 6, 7, 8, 9, 10, 11, // edges 365f34f2dcSJed Brown 0, 1, 3, 2, // centers 375f34f2dcSJed Brown }; 385f34f2dcSJed Brown static const int quad_25[] = { 395f34f2dcSJed Brown 21, 22, 23, 24, // vertices 405f34f2dcSJed Brown 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges 415f34f2dcSJed Brown 0, 1, 2, 5, 8, 7, 6, 3, 4, // centers 425f34f2dcSJed Brown }; 435f34f2dcSJed Brown static const int tetra_4[4] = {0, 2, 1, 3}; 445f34f2dcSJed Brown static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4}; 455f34f2dcSJed Brown static const int tetra_20[20] = { 465f34f2dcSJed Brown 16, 18, 17, 19, // vertices 475f34f2dcSJed Brown 9, 8, 7, 6, 5, 4, // bottom edges 485f34f2dcSJed Brown 10, 11, 14, 15, 13, 12, // side edges 495f34f2dcSJed Brown 0, 2, 3, 1, // faces 505f34f2dcSJed Brown }; 515f34f2dcSJed Brown static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 525f34f2dcSJed Brown static const int hexa_27[27] = { 535f34f2dcSJed Brown 19, 22, 21, 20, 23, 24, 25, 26, // vertices 545f34f2dcSJed Brown 10, 9, 8, 7, // bottom edges 555f34f2dcSJed Brown 16, 15, 18, 17, // mid edges 565f34f2dcSJed Brown 11, 12, 13, 14, // top edges 575f34f2dcSJed Brown 1, 3, 5, 4, 6, 2, // faces 585f34f2dcSJed Brown 0, // center 595f34f2dcSJed Brown }; 609371c9d4SSatish Balay static const int hexa_64[64] = { 619371c9d4SSatish Balay // 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 625f34f2dcSJed Brown 56, 59, 58, 57, 60, 61, 62, 63, // vertices 635f34f2dcSJed Brown 39, 38, 37, 36, 35, 34, 33, 32, // bottom edges 645f34f2dcSJed Brown 51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24 655f34f2dcSJed Brown 40, 41, 42, 43, 44, 45, 46, 47, // top edges 665f34f2dcSJed Brown 8, 10, 11, 9, // z-minus face 675f34f2dcSJed Brown 16, 17, 19, 18, // y-minus face 685f34f2dcSJed Brown 24, 25, 27, 26, // x-plus face 695f34f2dcSJed Brown 20, 21, 23, 22, // y-plus face 705f34f2dcSJed Brown 30, 28, 29, 31, // x-minus face 715f34f2dcSJed Brown 12, 13, 15, 14, // z-plus face 725f34f2dcSJed Brown 0, 1, 3, 2, 4, 5, 7, 6, // center 735f34f2dcSJed Brown }; 745f34f2dcSJed Brown 755f34f2dcSJed Brown PetscFunctionBegin; 76472b9844SJames Wright element_type_tmp = CGNS_ENUMV(ElementTypeNull); 775f34f2dcSJed Brown *perm = NULL; 785f34f2dcSJed Brown switch (cell_type) { 795f34f2dcSJed Brown case DM_POLYTOPE_SEGMENT: 805f34f2dcSJed Brown switch (closure_size) { 81d71ae5a4SJacob Faibussowitsch case 2: 82472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_2); 83d71ae5a4SJacob Faibussowitsch *perm = bar_2; 84d71ae5a4SJacob Faibussowitsch case 3: 85472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_3); 86d71ae5a4SJacob Faibussowitsch *perm = bar_3; 875f34f2dcSJed Brown case 4: 88472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_4); 895f34f2dcSJed Brown *perm = bar_4; 905f34f2dcSJed Brown break; 915f34f2dcSJed Brown case 5: 92472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_5); 935f34f2dcSJed Brown *perm = bar_5; 945f34f2dcSJed Brown break; 95d71ae5a4SJacob Faibussowitsch default: 96d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 975f34f2dcSJed Brown } 985f34f2dcSJed Brown break; 995f34f2dcSJed Brown case DM_POLYTOPE_TRIANGLE: 1005f34f2dcSJed Brown switch (closure_size) { 1015f34f2dcSJed Brown case 3: 102472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_3); 1035f34f2dcSJed Brown *perm = tri_3; 1045f34f2dcSJed Brown break; 1055f34f2dcSJed Brown case 6: 106472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_6); 1075f34f2dcSJed Brown *perm = tri_6; 1085f34f2dcSJed Brown break; 1095f34f2dcSJed Brown case 10: 110472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_10); 1115f34f2dcSJed Brown *perm = tri_10; 1125f34f2dcSJed Brown break; 113d71ae5a4SJacob Faibussowitsch default: 114d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1155f34f2dcSJed Brown } 1165f34f2dcSJed Brown break; 1175f34f2dcSJed Brown case DM_POLYTOPE_QUADRILATERAL: 1185f34f2dcSJed Brown switch (closure_size) { 1195f34f2dcSJed Brown case 4: 120472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_4); 1215f34f2dcSJed Brown *perm = quad_4; 1225f34f2dcSJed Brown break; 1235f34f2dcSJed Brown case 9: 124472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_9); 1255f34f2dcSJed Brown *perm = quad_9; 1265f34f2dcSJed Brown break; 1275f34f2dcSJed Brown case 16: 128472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_16); 1295f34f2dcSJed Brown *perm = quad_16; 1305f34f2dcSJed Brown break; 1315f34f2dcSJed Brown case 25: 132472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_25); 1335f34f2dcSJed Brown *perm = quad_25; 1345f34f2dcSJed Brown break; 135d71ae5a4SJacob Faibussowitsch default: 136d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1375f34f2dcSJed Brown } 1385f34f2dcSJed Brown break; 1395f34f2dcSJed Brown case DM_POLYTOPE_TETRAHEDRON: 1405f34f2dcSJed Brown switch (closure_size) { 1415f34f2dcSJed Brown case 4: 142472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_4); 1435f34f2dcSJed Brown *perm = tetra_4; 1445f34f2dcSJed Brown break; 1455f34f2dcSJed Brown case 10: 146472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_10); 1475f34f2dcSJed Brown *perm = tetra_10; 1485f34f2dcSJed Brown break; 1495f34f2dcSJed Brown case 20: 150472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_20); 1515f34f2dcSJed Brown *perm = tetra_20; 1525f34f2dcSJed Brown break; 153d71ae5a4SJacob Faibussowitsch default: 154d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1555f34f2dcSJed Brown } 1565f34f2dcSJed Brown break; 1575f34f2dcSJed Brown case DM_POLYTOPE_HEXAHEDRON: 1585f34f2dcSJed Brown switch (closure_size) { 1595f34f2dcSJed Brown case 8: 160472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_8); 1615f34f2dcSJed Brown *perm = hexa_8; 1625f34f2dcSJed Brown break; 1635f34f2dcSJed Brown case 27: 164472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_27); 1655f34f2dcSJed Brown *perm = hexa_27; 1665f34f2dcSJed Brown break; 1675f34f2dcSJed Brown case 64: 168472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_64); 1695f34f2dcSJed Brown *perm = hexa_64; 1705f34f2dcSJed Brown break; 171d71ae5a4SJacob Faibussowitsch default: 172d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1735f34f2dcSJed Brown } 1745f34f2dcSJed Brown break; 175d71ae5a4SJacob Faibussowitsch default: 176d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1775f34f2dcSJed Brown } 178472b9844SJames Wright if (element_type) *element_type = element_type_tmp; 179472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 180472b9844SJames Wright } 181472b9844SJames Wright 182472b9844SJames Wright /* 183472b9844SJames Wright Input Parameters: 184472b9844SJames Wright + cellType - The CGNS-defined element type 185472b9844SJames Wright 186472b9844SJames Wright Output Parameters: 187472b9844SJames Wright + dmcelltype - The equivalent DMPolytopeType for the cellType 188472b9844SJames Wright . numCorners - Number of corners of the polytope 189472b9844SJames Wright - dim - The topological dimension of the polytope 190472b9844SJames Wright 191472b9844SJames Wright CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 192472b9844SJames Wright */ 193472b9844SJames Wright static inline PetscErrorCode CGNSElementTypeGetTopologyInfo(CGNS_ENUMT(ElementType_t) cellType, DMPolytopeType *dmcelltype, PetscInt *numCorners, PetscInt *dim) 194472b9844SJames Wright { 195472b9844SJames Wright DMPolytopeType _dmcelltype; 196472b9844SJames Wright 197472b9844SJames Wright PetscFunctionBeginUser; 198472b9844SJames Wright switch (cellType) { 199472b9844SJames Wright case CGNS_ENUMV(BAR_2): 200472b9844SJames Wright case CGNS_ENUMV(BAR_3): 201472b9844SJames Wright case CGNS_ENUMV(BAR_4): 202472b9844SJames Wright case CGNS_ENUMV(BAR_5): 203472b9844SJames Wright _dmcelltype = DM_POLYTOPE_SEGMENT; 204472b9844SJames Wright break; 205472b9844SJames Wright case CGNS_ENUMV(TRI_3): 206472b9844SJames Wright case CGNS_ENUMV(TRI_6): 207472b9844SJames Wright case CGNS_ENUMV(TRI_9): 208472b9844SJames Wright case CGNS_ENUMV(TRI_10): 209472b9844SJames Wright case CGNS_ENUMV(TRI_12): 210472b9844SJames Wright case CGNS_ENUMV(TRI_15): 211472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TRIANGLE; 212472b9844SJames Wright break; 213472b9844SJames Wright case CGNS_ENUMV(QUAD_4): 214472b9844SJames Wright case CGNS_ENUMV(QUAD_8): 215472b9844SJames Wright case CGNS_ENUMV(QUAD_9): 216472b9844SJames Wright case CGNS_ENUMV(QUAD_12): 217472b9844SJames Wright case CGNS_ENUMV(QUAD_16): 218472b9844SJames Wright case CGNS_ENUMV(QUAD_P4_16): 219472b9844SJames Wright case CGNS_ENUMV(QUAD_25): 220472b9844SJames Wright _dmcelltype = DM_POLYTOPE_QUADRILATERAL; 221472b9844SJames Wright break; 222472b9844SJames Wright case CGNS_ENUMV(TETRA_4): 223472b9844SJames Wright case CGNS_ENUMV(TETRA_10): 224472b9844SJames Wright case CGNS_ENUMV(TETRA_16): 225472b9844SJames Wright case CGNS_ENUMV(TETRA_20): 226472b9844SJames Wright case CGNS_ENUMV(TETRA_22): 227472b9844SJames Wright case CGNS_ENUMV(TETRA_34): 228472b9844SJames Wright case CGNS_ENUMV(TETRA_35): 229472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TETRAHEDRON; 230472b9844SJames Wright break; 231472b9844SJames Wright case CGNS_ENUMV(PYRA_5): 232472b9844SJames Wright case CGNS_ENUMV(PYRA_13): 233472b9844SJames Wright case CGNS_ENUMV(PYRA_14): 234472b9844SJames Wright case CGNS_ENUMV(PYRA_21): 235472b9844SJames Wright case CGNS_ENUMV(PYRA_29): 236472b9844SJames Wright case CGNS_ENUMV(PYRA_P4_29): 237472b9844SJames Wright case CGNS_ENUMV(PYRA_30): 238472b9844SJames Wright case CGNS_ENUMV(PYRA_50): 239472b9844SJames Wright case CGNS_ENUMV(PYRA_55): 240472b9844SJames Wright _dmcelltype = DM_POLYTOPE_PYRAMID; 241472b9844SJames Wright break; 242472b9844SJames Wright case CGNS_ENUMV(PENTA_6): 243472b9844SJames Wright case CGNS_ENUMV(PENTA_15): 244472b9844SJames Wright case CGNS_ENUMV(PENTA_18): 245472b9844SJames Wright case CGNS_ENUMV(PENTA_24): 246472b9844SJames Wright case CGNS_ENUMV(PENTA_33): 247472b9844SJames Wright case CGNS_ENUMV(PENTA_38): 248472b9844SJames Wright case CGNS_ENUMV(PENTA_40): 249472b9844SJames Wright case CGNS_ENUMV(PENTA_66): 250472b9844SJames Wright case CGNS_ENUMV(PENTA_75): 251472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TRI_PRISM; 252472b9844SJames Wright break; 253472b9844SJames Wright case CGNS_ENUMV(HEXA_8): 254472b9844SJames Wright case CGNS_ENUMV(HEXA_20): 255472b9844SJames Wright case CGNS_ENUMV(HEXA_27): 256472b9844SJames Wright case CGNS_ENUMV(HEXA_32): 257472b9844SJames Wright case CGNS_ENUMV(HEXA_44): 258472b9844SJames Wright case CGNS_ENUMV(HEXA_56): 259472b9844SJames Wright case CGNS_ENUMV(HEXA_64): 260472b9844SJames Wright case CGNS_ENUMV(HEXA_98): 261472b9844SJames Wright case CGNS_ENUMV(HEXA_125): 262472b9844SJames Wright _dmcelltype = DM_POLYTOPE_HEXAHEDRON; 263472b9844SJames Wright break; 264472b9844SJames Wright case CGNS_ENUMV(MIXED): 265472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 266472b9844SJames Wright default: 267472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: %d", (int)cellType); 268472b9844SJames Wright } 269472b9844SJames Wright 270472b9844SJames Wright if (dmcelltype) *dmcelltype = _dmcelltype; 271472b9844SJames Wright if (numCorners) *numCorners = DMPolytopeTypeGetNumVertices(_dmcelltype); 272472b9844SJames Wright if (dim) *dim = DMPolytopeTypeGetDim(_dmcelltype); 273472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 274472b9844SJames Wright } 275472b9844SJames Wright 276472b9844SJames Wright /* 277472b9844SJames Wright Input Parameters: 278472b9844SJames Wright + cellType - The CGNS-defined cell type 279472b9844SJames Wright 280472b9844SJames Wright Output Parameters: 281472b9844SJames Wright + numClosure - Number of nodes that define the function space on the cell 282472b9844SJames Wright - pOrder - The polynomial order of the cell 283472b9844SJames Wright 284472b9844SJames Wright CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 285472b9844SJames Wright 286472b9844SJames Wright Note: we only support "full" elements, ie. not seredipity elements 287472b9844SJames Wright */ 288472b9844SJames Wright static inline PetscErrorCode CGNSElementTypeGetDiscretizationInfo(CGNS_ENUMT(ElementType_t) cellType, PetscInt *numClosure, PetscInt *pOrder) 289472b9844SJames Wright { 290472b9844SJames Wright PetscInt _numClosure, _pOrder; 291472b9844SJames Wright 292472b9844SJames Wright PetscFunctionBeginUser; 293472b9844SJames Wright switch (cellType) { 294472b9844SJames Wright case CGNS_ENUMV(BAR_2): 295472b9844SJames Wright _numClosure = 2; 296472b9844SJames Wright _pOrder = 1; 297472b9844SJames Wright break; 298472b9844SJames Wright case CGNS_ENUMV(BAR_3): 299472b9844SJames Wright _numClosure = 3; 300472b9844SJames Wright _pOrder = 2; 301472b9844SJames Wright break; 302472b9844SJames Wright case CGNS_ENUMV(BAR_4): 303472b9844SJames Wright _numClosure = 4; 304472b9844SJames Wright _pOrder = 3; 305472b9844SJames Wright break; 306472b9844SJames Wright case CGNS_ENUMV(BAR_5): 307472b9844SJames Wright _numClosure = 5; 308472b9844SJames Wright _pOrder = 4; 309472b9844SJames Wright break; 310472b9844SJames Wright case CGNS_ENUMV(TRI_3): 311472b9844SJames Wright _numClosure = 3; 312472b9844SJames Wright _pOrder = 1; 313472b9844SJames Wright break; 314472b9844SJames Wright case CGNS_ENUMV(TRI_6): 315472b9844SJames Wright _numClosure = 6; 316472b9844SJames Wright _pOrder = 2; 317472b9844SJames Wright break; 318472b9844SJames Wright case CGNS_ENUMV(TRI_10): 319472b9844SJames Wright _numClosure = 10; 320472b9844SJames Wright _pOrder = 3; 321472b9844SJames Wright break; 322472b9844SJames Wright case CGNS_ENUMV(TRI_15): 323472b9844SJames Wright _numClosure = 15; 324472b9844SJames Wright _pOrder = 4; 325472b9844SJames Wright break; 326472b9844SJames Wright case CGNS_ENUMV(QUAD_4): 327472b9844SJames Wright _numClosure = 4; 328472b9844SJames Wright _pOrder = 1; 329472b9844SJames Wright break; 330472b9844SJames Wright case CGNS_ENUMV(QUAD_9): 331472b9844SJames Wright _numClosure = 9; 332472b9844SJames Wright _pOrder = 2; 333472b9844SJames Wright break; 334472b9844SJames Wright case CGNS_ENUMV(QUAD_16): 335472b9844SJames Wright _numClosure = 16; 336472b9844SJames Wright _pOrder = 3; 337472b9844SJames Wright break; 338472b9844SJames Wright case CGNS_ENUMV(QUAD_25): 339472b9844SJames Wright _numClosure = 25; 340472b9844SJames Wright _pOrder = 4; 341472b9844SJames Wright break; 342472b9844SJames Wright case CGNS_ENUMV(TETRA_4): 343472b9844SJames Wright _numClosure = 4; 344472b9844SJames Wright _pOrder = 1; 345472b9844SJames Wright break; 346472b9844SJames Wright case CGNS_ENUMV(TETRA_10): 347472b9844SJames Wright _numClosure = 10; 348472b9844SJames Wright _pOrder = 2; 349472b9844SJames Wright break; 350472b9844SJames Wright case CGNS_ENUMV(TETRA_20): 351472b9844SJames Wright _numClosure = 20; 352472b9844SJames Wright _pOrder = 3; 353472b9844SJames Wright break; 354472b9844SJames Wright case CGNS_ENUMV(TETRA_35): 355472b9844SJames Wright _numClosure = 35; 356472b9844SJames Wright _pOrder = 4; 357472b9844SJames Wright break; 358472b9844SJames Wright case CGNS_ENUMV(PYRA_5): 359472b9844SJames Wright _numClosure = 5; 360472b9844SJames Wright _pOrder = 1; 361472b9844SJames Wright break; 362472b9844SJames Wright case CGNS_ENUMV(PYRA_14): 363472b9844SJames Wright _numClosure = 14; 364472b9844SJames Wright _pOrder = 2; 365472b9844SJames Wright break; 366472b9844SJames Wright case CGNS_ENUMV(PYRA_30): 367472b9844SJames Wright _numClosure = 30; 368472b9844SJames Wright _pOrder = 3; 369472b9844SJames Wright break; 370472b9844SJames Wright case CGNS_ENUMV(PYRA_55): 371472b9844SJames Wright _numClosure = 55; 372472b9844SJames Wright _pOrder = 4; 373472b9844SJames Wright break; 374472b9844SJames Wright case CGNS_ENUMV(PENTA_6): 375472b9844SJames Wright _numClosure = 6; 376472b9844SJames Wright _pOrder = 1; 377472b9844SJames Wright break; 378472b9844SJames Wright case CGNS_ENUMV(PENTA_18): 379472b9844SJames Wright _numClosure = 18; 380472b9844SJames Wright _pOrder = 2; 381472b9844SJames Wright break; 382472b9844SJames Wright case CGNS_ENUMV(PENTA_40): 383472b9844SJames Wright _numClosure = 40; 384472b9844SJames Wright _pOrder = 3; 385472b9844SJames Wright break; 386472b9844SJames Wright case CGNS_ENUMV(PENTA_75): 387472b9844SJames Wright _numClosure = 75; 388472b9844SJames Wright _pOrder = 4; 389472b9844SJames Wright break; 390472b9844SJames Wright case CGNS_ENUMV(HEXA_8): 391472b9844SJames Wright _numClosure = 8; 392472b9844SJames Wright _pOrder = 1; 393472b9844SJames Wright break; 394472b9844SJames Wright case CGNS_ENUMV(HEXA_27): 395472b9844SJames Wright _numClosure = 27; 396472b9844SJames Wright _pOrder = 2; 397472b9844SJames Wright break; 398472b9844SJames Wright case CGNS_ENUMV(HEXA_64): 399472b9844SJames Wright _numClosure = 64; 400472b9844SJames Wright _pOrder = 3; 401472b9844SJames Wright break; 402472b9844SJames Wright case CGNS_ENUMV(HEXA_125): 403472b9844SJames Wright _numClosure = 125; 404472b9844SJames Wright _pOrder = 4; 405472b9844SJames Wright break; 406472b9844SJames Wright case CGNS_ENUMV(MIXED): 407472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 408472b9844SJames Wright default: 409472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported or Invalid cell type %d", (int)cellType); 410472b9844SJames Wright } 411472b9844SJames Wright if (numClosure) *numClosure = _numClosure; 412472b9844SJames Wright if (pOrder) *pOrder = _pOrder; 413472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 414472b9844SJames Wright } 415472b9844SJames Wright 416472b9844SJames Wright static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd) 417472b9844SJames Wright { 418472b9844SJames Wright PetscFunctionBegin; 419472b9844SJames Wright switch (pd) { 420472b9844SJames Wright case PETSC_FLOAT: 421472b9844SJames Wright *cd = CGNS_ENUMV(RealSingle); 422472b9844SJames Wright break; 423472b9844SJames Wright case PETSC_DOUBLE: 424472b9844SJames Wright *cd = CGNS_ENUMV(RealDouble); 425472b9844SJames Wright break; 426472b9844SJames Wright case PETSC_COMPLEX: 427472b9844SJames Wright *cd = CGNS_ENUMV(ComplexDouble); 428472b9844SJames Wright break; 429472b9844SJames Wright default: 430472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]); 431472b9844SJames Wright } 432472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 433472b9844SJames Wright } 434472b9844SJames Wright 435472b9844SJames Wright PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 436472b9844SJames Wright { 437472b9844SJames Wright int cgid = -1; 438472b9844SJames Wright PetscBool use_parallel_viewer = PETSC_FALSE; 439472b9844SJames Wright 440472b9844SJames Wright PetscFunctionBegin; 441472b9844SJames Wright PetscAssertPointer(filename, 2); 4425582b114SJames Wright PetscCall(PetscViewerCGNSRegisterLogEvents_Internal()); 443472b9844SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 444472b9844SJames Wright 445472b9844SJames Wright if (use_parallel_viewer) { 446472b9844SJames Wright PetscCallCGNS(cgp_mpi_comm(comm)); 4474c155f28SJames Wright PetscCallCGNSOpen(cgp_open(filename, CG_MODE_READ, &cgid), 0, 0); 448472b9844SJames Wright PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cgp_open(\"%s\",...) did not return a valid file ID", filename); 449472b9844SJames Wright PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 4504c155f28SJames Wright PetscCallCGNSClose(cgp_close(cgid), 0, 0); 451472b9844SJames Wright } else { 4524c155f28SJames Wright PetscCallCGNSOpen(cg_open(filename, CG_MODE_READ, &cgid), 0, 0); 453472b9844SJames Wright PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 454472b9844SJames Wright PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 4554c155f28SJames Wright PetscCallCGNSClose(cg_close(cgid), 0, 0); 456472b9844SJames Wright } 457472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 458472b9844SJames Wright } 459472b9844SJames Wright 460472b9844SJames Wright PetscErrorCode DMPlexCreateCGNS_Internal_Serial(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 461472b9844SJames Wright { 462472b9844SJames Wright PetscMPIInt num_proc, rank; 463472b9844SJames Wright DM cdm; 464472b9844SJames Wright DMLabel label; 465472b9844SJames Wright PetscSection coordSection; 466472b9844SJames Wright Vec coordinates; 467472b9844SJames Wright PetscScalar *coords; 468472b9844SJames Wright PetscInt *cellStart, *vertStart, v; 469472b9844SJames Wright PetscInt labelIdRange[2], labelId; 470472b9844SJames Wright /* Read from file */ 471472b9844SJames Wright char basename[CGIO_MAX_NAME_LENGTH + 1]; 472472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 473472b9844SJames Wright int dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0; 474472b9844SJames Wright int nzones = 0; 475472b9844SJames Wright const int B = 1; // Only support single base 476472b9844SJames Wright 477472b9844SJames Wright PetscFunctionBegin; 478472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 479472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 480472b9844SJames Wright PetscCall(DMCreate(comm, dm)); 481472b9844SJames Wright PetscCall(DMSetType(*dm, DMPLEX)); 482472b9844SJames Wright 483472b9844SJames Wright /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 484472b9844SJames Wright if (rank == 0) { 485472b9844SJames Wright int nbases, z; 486472b9844SJames Wright 4874c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 488472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 4894c155f28SJames Wright PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0); 4904c155f28SJames Wright PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0); 491472b9844SJames Wright PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart)); 492472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 493472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 494472b9844SJames Wright 4954c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 496472b9844SJames Wright numVertices += sizes[0]; 497472b9844SJames Wright numCells += sizes[1]; 498472b9844SJames Wright cellStart[z] += sizes[1] + cellStart[z - 1]; 499472b9844SJames Wright vertStart[z] += sizes[0] + vertStart[z - 1]; 500472b9844SJames Wright } 501472b9844SJames Wright for (z = 1; z <= nzones; ++z) vertStart[z] += numCells; 502472b9844SJames Wright coordDim = dim; 503472b9844SJames Wright } 504472b9844SJames Wright PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm)); 505472b9844SJames Wright PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 506472b9844SJames Wright PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm)); 507472b9844SJames Wright PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm)); 508472b9844SJames Wright 509472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 510472b9844SJames Wright PetscCall(DMSetDimension(*dm, dim)); 511472b9844SJames Wright PetscCall(DMCreateLabel(*dm, "celltype")); 512472b9844SJames Wright PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 513472b9844SJames Wright 514472b9844SJames Wright /* Read zone information */ 515472b9844SJames Wright if (rank == 0) { 516472b9844SJames Wright int z, c, c_loc; 517472b9844SJames Wright 518472b9844SJames Wright /* Read the cell set connectivity table and build mesh topology 519472b9844SJames Wright CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 520472b9844SJames Wright /* First set sizes */ 521472b9844SJames Wright for (z = 1, c = 0; z <= nzones; ++z) { 522472b9844SJames Wright CGNS_ENUMT(ZoneType_t) zonetype; 523472b9844SJames Wright int nsections; 524472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 525472b9844SJames Wright cgsize_t start, end; 526472b9844SJames Wright int nbndry, parentFlag; 527472b9844SJames Wright PetscInt numCorners, pOrder; 528472b9844SJames Wright DMPolytopeType ctype; 529472b9844SJames Wright const int S = 1; // Only support single section 530472b9844SJames Wright 5314c155f28SJames Wright PetscCallCGNSRead(cg_zone_type(cgid, B, z, &zonetype), *dm, 0); 532472b9844SJames Wright PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS"); 5334c155f28SJames Wright PetscCallCGNSRead(cg_nsections(cgid, B, z, &nsections), *dm, 0); 534472b9844SJames Wright PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections); 5354c155f28SJames Wright PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 536472b9844SJames Wright if (cellType == CGNS_ENUMV(MIXED)) { 537472b9844SJames Wright cgsize_t elementDataSize, *elements; 538472b9844SJames Wright PetscInt off; 539472b9844SJames Wright 5404c155f28SJames Wright PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 541472b9844SJames Wright PetscCall(PetscMalloc1(elementDataSize, &elements)); 5424c155f28SJames Wright PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 543472b9844SJames Wright for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 544ca328833SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[off], &ctype, &numCorners, NULL)); 545ca328833SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[off], NULL, &pOrder)); 546472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 547472b9844SJames Wright PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 548472b9844SJames Wright PetscCall(DMPlexSetCellType(*dm, c, ctype)); 549472b9844SJames Wright off += numCorners + 1; 550472b9844SJames Wright } 551472b9844SJames Wright PetscCall(PetscFree(elements)); 552472b9844SJames Wright } else { 553472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &ctype, &numCorners, NULL)); 554472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 555472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 556472b9844SJames Wright for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 557472b9844SJames Wright PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 558472b9844SJames Wright PetscCall(DMPlexSetCellType(*dm, c, ctype)); 559472b9844SJames Wright } 560472b9844SJames Wright } 561472b9844SJames Wright } 562472b9844SJames Wright for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 563472b9844SJames Wright } 564472b9844SJames Wright 565472b9844SJames Wright PetscCall(DMSetUp(*dm)); 566472b9844SJames Wright 567472b9844SJames Wright PetscCall(DMCreateLabel(*dm, "zone")); 568472b9844SJames Wright if (rank == 0) { 569472b9844SJames Wright int z, c, c_loc, v_loc; 570472b9844SJames Wright 571472b9844SJames Wright PetscCall(DMGetLabel(*dm, "zone", &label)); 572472b9844SJames Wright for (z = 1, c = 0; z <= nzones; ++z) { 573472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 574472b9844SJames Wright cgsize_t elementDataSize, *elements, start, end; 575472b9844SJames Wright int nbndry, parentFlag; 576472b9844SJames Wright PetscInt *cone, numc, numCorners, maxCorners = 27, pOrder; 577472b9844SJames Wright const int S = 1; // Only support single section 578472b9844SJames Wright 5794c155f28SJames Wright PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 580472b9844SJames Wright numc = end - start; 5814c155f28SJames Wright PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 582472b9844SJames Wright PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone)); 5834c155f28SJames Wright PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 584472b9844SJames Wright if (cellType == CGNS_ENUMV(MIXED)) { 585472b9844SJames Wright /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 586472b9844SJames Wright for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 587ca328833SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &numCorners, NULL)); 588ca328833SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &pOrder)); 589472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 590472b9844SJames Wright ++v; 591472b9844SJames Wright for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 592472b9844SJames Wright PetscCall(DMPlexReorderCell(*dm, c, cone)); 593472b9844SJames Wright PetscCall(DMPlexSetCone(*dm, c, cone)); 594472b9844SJames Wright PetscCall(DMLabelSetValue(label, c, z)); 595472b9844SJames Wright } 596472b9844SJames Wright } else { 597472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, NULL, &numCorners, NULL)); 598472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 599472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 600472b9844SJames Wright /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 601472b9844SJames Wright for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 602472b9844SJames Wright for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 603472b9844SJames Wright PetscCall(DMPlexReorderCell(*dm, c, cone)); 604472b9844SJames Wright PetscCall(DMPlexSetCone(*dm, c, cone)); 605472b9844SJames Wright PetscCall(DMLabelSetValue(label, c, z)); 606472b9844SJames Wright } 607472b9844SJames Wright } 608472b9844SJames Wright PetscCall(PetscFree2(elements, cone)); 609472b9844SJames Wright } 610472b9844SJames Wright } 611472b9844SJames Wright 612472b9844SJames Wright PetscCall(DMPlexSymmetrize(*dm)); 613472b9844SJames Wright PetscCall(DMPlexStratify(*dm)); 61401432a30SJames Wright if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 615472b9844SJames Wright 616472b9844SJames Wright /* Read coordinates */ 617472b9844SJames Wright PetscCall(DMSetCoordinateDim(*dm, coordDim)); 618472b9844SJames Wright PetscCall(DMGetCoordinateDM(*dm, &cdm)); 619472b9844SJames Wright PetscCall(DMGetLocalSection(cdm, &coordSection)); 620472b9844SJames Wright PetscCall(PetscSectionSetNumFields(coordSection, 1)); 621472b9844SJames Wright PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim)); 622472b9844SJames Wright PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 623472b9844SJames Wright for (v = numCells; v < numCells + numVertices; ++v) { 624472b9844SJames Wright PetscCall(PetscSectionSetDof(coordSection, v, dim)); 625472b9844SJames Wright PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim)); 626472b9844SJames Wright } 627472b9844SJames Wright PetscCall(PetscSectionSetUp(coordSection)); 628472b9844SJames Wright 629472b9844SJames Wright PetscCall(DMCreateLocalVector(cdm, &coordinates)); 630472b9844SJames Wright PetscCall(VecGetArray(coordinates, &coords)); 631472b9844SJames Wright if (rank == 0) { 632472b9844SJames Wright PetscInt off = 0; 633472b9844SJames Wright float *x[3]; 634472b9844SJames Wright int z, d; 635472b9844SJames Wright 636472b9844SJames Wright PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2])); 637472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 638472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 639472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 640472b9844SJames Wright cgsize_t range_min[3] = {1, 1, 1}; 641472b9844SJames Wright cgsize_t range_max[3] = {1, 1, 1}; 642472b9844SJames Wright int ngrids, ncoords; 643472b9844SJames Wright 6444c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 645472b9844SJames Wright range_max[0] = sizes[0]; 6464c155f28SJames Wright PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0); 647472b9844SJames Wright PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 6484c155f28SJames Wright PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0); 649472b9844SJames Wright PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 650472b9844SJames Wright for (d = 0; d < coordDim; ++d) { 6514c155f28SJames Wright PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1 + d, &datatype, buffer), *dm, 0); 6524c155f28SJames Wright PetscCallCGNSReadData(cg_coord_read(cgid, B, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]), *dm, 0); 653472b9844SJames Wright } 654472b9844SJames Wright if (coordDim >= 1) { 655472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v]; 656472b9844SJames Wright } 657472b9844SJames Wright if (coordDim >= 2) { 658472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v]; 659472b9844SJames Wright } 660472b9844SJames Wright if (coordDim >= 3) { 661472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v]; 662472b9844SJames Wright } 663472b9844SJames Wright off += sizes[0]; 664472b9844SJames Wright } 665472b9844SJames Wright PetscCall(PetscFree3(x[0], x[1], x[2])); 666472b9844SJames Wright } 667472b9844SJames Wright PetscCall(VecRestoreArray(coordinates, &coords)); 668472b9844SJames Wright 669472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 670472b9844SJames Wright PetscCall(VecSetBlockSize(coordinates, coordDim)); 671472b9844SJames Wright PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 672472b9844SJames Wright PetscCall(VecDestroy(&coordinates)); 673472b9844SJames Wright 674472b9844SJames Wright /* Read boundary conditions */ 675472b9844SJames Wright PetscCall(DMGetNumLabels(*dm, &labelIdRange[0])); 676472b9844SJames Wright if (rank == 0) { 677472b9844SJames Wright CGNS_ENUMT(BCType_t) bctype; 678472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 679472b9844SJames Wright CGNS_ENUMT(PointSetType_t) pointtype; 680472b9844SJames Wright cgsize_t *points; 681472b9844SJames Wright PetscReal *normals; 682472b9844SJames Wright int normal[3]; 683472b9844SJames Wright char *bcname = buffer; 684472b9844SJames Wright cgsize_t npoints, nnormals; 685472b9844SJames Wright int z, nbc, bc, c, ndatasets; 686472b9844SJames Wright 687472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 6884c155f28SJames Wright PetscCallCGNSRead(cg_nbocos(cgid, B, z, &nbc), *dm, 0); 689472b9844SJames Wright for (bc = 1; bc <= nbc; ++bc) { 6904c155f28SJames Wright PetscCallCGNSRead(cg_boco_info(cgid, B, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets), *dm, 0); 691472b9844SJames Wright PetscCall(DMCreateLabel(*dm, bcname)); 692472b9844SJames Wright PetscCall(DMGetLabel(*dm, bcname, &label)); 693472b9844SJames Wright PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals)); 6944c155f28SJames Wright PetscCallCGNSReadData(cg_boco_read(cgid, B, z, bc, points, (void *)normals), *dm, 0); 695472b9844SJames Wright if (pointtype == CGNS_ENUMV(ElementRange)) { 696472b9844SJames Wright // Range of cells: assuming half-open interval 697472b9844SJames Wright for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 698472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(ElementList)) { 699472b9844SJames Wright // List of cells 700472b9844SJames Wright for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 701472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(PointRange)) { 702472b9844SJames Wright CGNS_ENUMT(GridLocation_t) gridloc; 703472b9844SJames Wright 704472b9844SJames Wright // List of points: 705472b9844SJames Wright PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 7064c155f28SJames Wright PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 707472b9844SJames Wright // Range of points: assuming half-open interval 708472b9844SJames Wright for (c = points[0]; c < points[1]; ++c) { 709472b9844SJames Wright if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1)); 710472b9844SJames Wright else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 711472b9844SJames Wright } 712472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(PointList)) { 713472b9844SJames Wright CGNS_ENUMT(GridLocation_t) gridloc; 714472b9844SJames Wright 715472b9844SJames Wright // List of points: 716472b9844SJames Wright PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 7174c155f28SJames Wright PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 718472b9844SJames Wright for (c = 0; c < npoints; ++c) { 719472b9844SJames Wright if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1)); 720472b9844SJames Wright else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 721472b9844SJames Wright } 722472b9844SJames Wright } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype); 723472b9844SJames Wright PetscCall(PetscFree2(points, normals)); 724472b9844SJames Wright } 725472b9844SJames Wright } 726472b9844SJames Wright PetscCall(PetscFree2(cellStart, vertStart)); 727472b9844SJames Wright } 728472b9844SJames Wright PetscCall(DMGetNumLabels(*dm, &labelIdRange[1])); 729472b9844SJames Wright PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm)); 730472b9844SJames Wright 731472b9844SJames Wright /* Create BC labels at all processes */ 732472b9844SJames Wright for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 733472b9844SJames Wright char *labelName = buffer; 734472b9844SJames Wright size_t len = sizeof(buffer); 735472b9844SJames Wright const char *locName; 736472b9844SJames Wright 737472b9844SJames Wright if (rank == 0) { 738472b9844SJames Wright PetscCall(DMGetLabelByNum(*dm, labelId, &label)); 739472b9844SJames Wright PetscCall(PetscObjectGetName((PetscObject)label, &locName)); 740472b9844SJames Wright PetscCall(PetscStrncpy(labelName, locName, len)); 741472b9844SJames Wright } 742472b9844SJames Wright PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm)); 743472b9844SJames Wright PetscCallMPI(DMCreateLabel(*dm, labelName)); 744472b9844SJames Wright } 745472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 746472b9844SJames Wright } 747472b9844SJames Wright 748472b9844SJames Wright PetscErrorCode DMPlexCreateCGNS_Internal_Parallel(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 749472b9844SJames Wright { 750472b9844SJames Wright PetscMPIInt num_proc, rank; 751472b9844SJames Wright /* Read from file */ 752472b9844SJames Wright char basename[CGIO_MAX_NAME_LENGTH + 1]; 753472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 754472b9844SJames Wright int dim = 0, physDim = 0, coordDim = 0; 755472b9844SJames Wright PetscInt NVertices = 0, NCells = 0; 756472b9844SJames Wright int nzones = 0, nbases; 757472b9844SJames Wright int z = 1; // Only supports single zone files 758472b9844SJames Wright int B = 1; // Only supports single base 759472b9844SJames Wright 760472b9844SJames Wright PetscFunctionBegin; 761472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 762472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 763472b9844SJames Wright PetscCall(DMCreate(comm, dm)); 764472b9844SJames Wright PetscCall(DMSetType(*dm, DMPLEX)); 765472b9844SJames Wright 7664c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 767472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 768472b9844SJames Wright // From the CGNS web page cell_dim phys_dim (embedding space in PETSC) CGNS defines as length of spatial vectors/components) 7694c155f28SJames Wright PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0); 7704c155f28SJames Wright PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0); 771472b9844SJames Wright PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Parallel reader limited to one zone, not %d", nzones); 772472b9844SJames Wright { 773472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 774472b9844SJames Wright 7754c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 776472b9844SJames Wright NVertices = sizes[0]; 777472b9844SJames Wright NCells = sizes[1]; 778472b9844SJames Wright } 779472b9844SJames Wright 780472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 781472b9844SJames Wright PetscCall(DMSetDimension(*dm, dim)); 782472b9844SJames Wright coordDim = dim; 783472b9844SJames Wright 784472b9844SJames Wright // 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 785472b9844SJames Wright // call to get this right but continuing for now with single section, single topology, one zone. 786472b9844SJames Wright // establish element ranges for my rank 787472b9844SJames Wright PetscInt mystarte, myende, mystartv, myendv, myownede, myownedv; 788472b9844SJames Wright PetscLayout elem_map, vtx_map; 789472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NCells, 1, &elem_map)); 790472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 791472b9844SJames Wright PetscCall(PetscLayoutGetRange(elem_map, &mystarte, &myende)); 792472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(elem_map, &myownede)); 793472b9844SJames Wright PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 794472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 795472b9844SJames Wright 796472b9844SJames Wright // -- Build Plex in parallel 797f1a88294SStefano Zampini DMPolytopeType dm_cell_type = DM_POLYTOPE_UNKNOWN; 798472b9844SJames Wright PetscInt pOrder = 1, numClosure = -1; 799472b9844SJames Wright cgsize_t *elements; 800472b9844SJames Wright { 801472b9844SJames Wright int nsections; 802472b9844SJames Wright PetscInt *elementsQ1, numCorners = -1; 803472b9844SJames Wright const int *perm; 804472b9844SJames Wright cgsize_t start, end; // Throwaway 805472b9844SJames Wright 806472b9844SJames Wright cg_nsections(cgid, B, z, &nsections); 807472b9844SJames Wright // Read element connectivity 808472b9844SJames Wright for (int index_sect = 1; index_sect <= nsections; index_sect++) { 809472b9844SJames Wright int nbndry, parentFlag; 810472b9844SJames Wright PetscInt cell_dim; 811472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 812472b9844SJames Wright 8134c155f28SJames Wright PetscCallCGNSRead(cg_section_read(cgid, B, z, index_sect, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 814472b9844SJames Wright 815472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &dm_cell_type, &numCorners, &cell_dim)); 816472b9844SJames Wright // Skip over element that are not max dimension (ie. boundary elements) 817472b9844SJames Wright if (cell_dim != dim) continue; 818472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, &numClosure, &pOrder)); 819472b9844SJames Wright PetscCall(PetscMalloc1(myownede * numClosure, &elements)); 8204c155f28SJames Wright PetscCallCGNSReadData(cgp_elements_read_data(cgid, B, z, index_sect, mystarte + 1, myende, elements), *dm, 0); 821472b9844SJames Wright for (PetscInt v = 0; v < myownede * numClosure; ++v) elements[v] -= 1; // 0 based 822472b9844SJames Wright break; 823472b9844SJames Wright } 824472b9844SJames Wright 825472b9844SJames Wright // Create corners-only connectivity 826472b9844SJames Wright PetscCall(PetscMalloc1(myownede * numCorners, &elementsQ1)); 827472b9844SJames Wright PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numCorners, NULL, &perm)); 828472b9844SJames Wright for (PetscInt e = 0; e < myownede; ++e) { 829472b9844SJames Wright for (PetscInt v = 0; v < numCorners; ++v) elementsQ1[e * numCorners + perm[v]] = elements[e * numClosure + v]; 830472b9844SJames Wright } 831472b9844SJames Wright 832472b9844SJames Wright // Build cell-vertex Plex 833472b9844SJames Wright PetscCall(DMPlexBuildFromCellListParallel(*dm, myownede, myownedv, NVertices, numCorners, elementsQ1, NULL, NULL)); 834472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-corner_dm_view")); 835472b9844SJames Wright PetscCall(PetscFree(elementsQ1)); 836472b9844SJames Wright } 837472b9844SJames Wright 83801432a30SJames Wright if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 839472b9844SJames Wright 840472b9844SJames Wright // -- Create SF for naive nodal-data read to elements 841472b9844SJames Wright PetscSF plex_to_cgns_sf; 842472b9844SJames Wright { 843472b9844SJames Wright PetscInt nleaves, num_comp; 844472b9844SJames Wright PetscInt *leaf, num_leaves = 0; 845472b9844SJames Wright PetscInt cStart, cEnd; 846472b9844SJames Wright const int *perm; 847472b9844SJames Wright PetscSF cgns_to_local_sf; 848472b9844SJames Wright PetscSection local_section; 849472b9844SJames Wright PetscFE fe; 850472b9844SJames Wright 851472b9844SJames Wright // sfNatural requires PetscSection to handle DMDistribute, so we use PetscFE to define the section 852472b9844SJames Wright // Use number of components = 1 to work with just the nodes themselves 853472b9844SJames Wright PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, dm_cell_type, pOrder, PETSC_DETERMINE, &fe)); 854472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "FE for sfNatural")); 855472b9844SJames Wright PetscCall(DMAddField(*dm, NULL, (PetscObject)fe)); 856472b9844SJames Wright PetscCall(DMCreateDS(*dm)); 857472b9844SJames Wright PetscCall(PetscFEDestroy(&fe)); 858472b9844SJames Wright 859472b9844SJames Wright PetscCall(DMGetLocalSection(*dm, &local_section)); 860472b9844SJames Wright PetscCall(PetscSectionViewFromOptions(local_section, NULL, "-fe_natural_section_view")); 861472b9844SJames Wright PetscCall(PetscSectionGetFieldComponents(local_section, 0, &num_comp)); 862472b9844SJames Wright PetscCall(PetscSectionGetStorageSize(local_section, &nleaves)); 863472b9844SJames Wright nleaves /= num_comp; 864472b9844SJames Wright PetscCall(PetscMalloc1(nleaves, &leaf)); 865472b9844SJames Wright for (PetscInt i = 0; i < nleaves; i++) leaf[i] = -1; 866472b9844SJames Wright 867472b9844SJames Wright // Get the permutation from CGNS closure numbering to PLEX closure numbering 868472b9844SJames Wright PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numClosure, NULL, &perm)); 869472b9844SJames Wright PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 870472b9844SJames Wright for (PetscInt cell = cStart; cell < cEnd; ++cell) { 871472b9844SJames Wright PetscInt num_closure_dof, *closure_idx = NULL; 872472b9844SJames Wright 873472b9844SJames Wright PetscCall(DMPlexGetClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 874472b9844SJames Wright PetscAssert(numClosure * num_comp == num_closure_dof, comm, PETSC_ERR_PLIB, "Closure dof size does not match polytope"); 875472b9844SJames Wright for (PetscInt i = 0; i < numClosure; i++) { 876472b9844SJames Wright PetscInt li = closure_idx[perm[i] * num_comp] / num_comp; 877472b9844SJames Wright if (li < 0) continue; 878472b9844SJames Wright 879472b9844SJames Wright PetscInt cgns_idx = elements[cell * numClosure + i]; 880472b9844SJames Wright if (leaf[li] == -1) { 881472b9844SJames Wright leaf[li] = cgns_idx; 882472b9844SJames Wright num_leaves++; 883472b9844SJames Wright } else PetscAssert(leaf[li] == cgns_idx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf does not match previously set"); 884472b9844SJames Wright } 885472b9844SJames Wright PetscCall(DMPlexRestoreClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 886472b9844SJames Wright } 887472b9844SJames Wright PetscAssert(num_leaves == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf count in closure does not match nleaves"); 888472b9844SJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)*dm), &cgns_to_local_sf)); 889472b9844SJames Wright PetscCall(PetscSFSetGraphLayout(cgns_to_local_sf, vtx_map, nleaves, NULL, PETSC_USE_POINTER, leaf)); 890472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)cgns_to_local_sf, "CGNS to Plex SF")); 891472b9844SJames Wright PetscCall(PetscSFViewFromOptions(cgns_to_local_sf, NULL, "-CGNStoPlex_sf_view")); 892472b9844SJames Wright PetscCall(PetscFree(leaf)); 893472b9844SJames Wright PetscCall(PetscFree(elements)); 894472b9844SJames Wright 895472b9844SJames Wright { // Convert cgns_to_local to global_to_cgns 896472b9844SJames Wright PetscSF sectionsf, cgns_to_global_sf; 897472b9844SJames Wright 898472b9844SJames Wright PetscCall(DMGetSectionSF(*dm, §ionsf)); 899472b9844SJames Wright PetscCall(PetscSFComposeInverse(cgns_to_local_sf, sectionsf, &cgns_to_global_sf)); 900472b9844SJames Wright PetscCall(PetscSFDestroy(&cgns_to_local_sf)); 901472b9844SJames Wright PetscCall(PetscSFCreateInverseSF(cgns_to_global_sf, &plex_to_cgns_sf)); 902472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)plex_to_cgns_sf, "Global Plex to CGNS")); 903472b9844SJames Wright PetscCall(PetscSFDestroy(&cgns_to_global_sf)); 904472b9844SJames Wright } 905472b9844SJames Wright } 906472b9844SJames Wright 907472b9844SJames Wright { // -- Set coordinates for DM 908472b9844SJames Wright PetscScalar *coords; 909472b9844SJames Wright float *x[3]; 910472b9844SJames Wright double *xd[3]; 911472b9844SJames Wright PetscBool read_with_double; 912472b9844SJames Wright PetscFE cfe; 913472b9844SJames Wright 914472b9844SJames Wright // Setup coordinate space first. Use pOrder here for isoparametric; revisit with CPEX-0045 High Order. 915472b9844SJames Wright PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, coordDim, dm_cell_type, pOrder, PETSC_DETERMINE, &cfe)); 9164c712d99Sksagiyam PetscCall(DMSetCoordinateDisc(*dm, cfe, PETSC_FALSE, PETSC_FALSE)); 917472b9844SJames Wright PetscCall(PetscFEDestroy(&cfe)); 918472b9844SJames Wright 919472b9844SJames Wright { // Determine if coords are written in single or double precision 920472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 921472b9844SJames Wright 9224c155f28SJames Wright PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1, &datatype, buffer), *dm, 0); 9235582b114SJames Wright read_with_double = datatype == CGNS_ENUMV(RealDouble) ? PETSC_TRUE : PETSC_FALSE; 924472b9844SJames Wright } 925472b9844SJames Wright 926472b9844SJames Wright // Read coords from file and set into component-major ordering 927472b9844SJames Wright if (read_with_double) PetscCall(PetscMalloc3(myownedv, &xd[0], myownedv, &xd[1], myownedv, &xd[2])); 928472b9844SJames Wright else PetscCall(PetscMalloc3(myownedv, &x[0], myownedv, &x[1], myownedv, &x[2])); 929472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * coordDim, &coords)); 930472b9844SJames Wright { 931472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 932472b9844SJames Wright cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 933472b9844SJames Wright cgsize_t range_max[3] = {myendv, 1, 1}; 934472b9844SJames Wright int ngrids, ncoords; 935472b9844SJames Wright 9364c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 9374c155f28SJames Wright PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0); 938472b9844SJames Wright PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 9394c155f28SJames Wright PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0); 940472b9844SJames Wright PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 941472b9844SJames Wright if (read_with_double) { 9424c155f28SJames Wright 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); 943472b9844SJames Wright if (coordDim >= 1) { 944472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = xd[0][v]; 945472b9844SJames Wright } 946472b9844SJames Wright if (coordDim >= 2) { 947472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = xd[1][v]; 948472b9844SJames Wright } 949472b9844SJames Wright if (coordDim >= 3) { 950472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = xd[2][v]; 951472b9844SJames Wright } 952472b9844SJames Wright } else { 9534c155f28SJames Wright 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); 954472b9844SJames Wright if (coordDim >= 1) { 955472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = x[0][v]; 956472b9844SJames Wright } 957472b9844SJames Wright if (coordDim >= 2) { 958472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = x[1][v]; 959472b9844SJames Wright } 960472b9844SJames Wright if (coordDim >= 3) { 961472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = x[2][v]; 962472b9844SJames Wright } 963472b9844SJames Wright } 964472b9844SJames Wright } 965472b9844SJames Wright if (read_with_double) PetscCall(PetscFree3(xd[0], xd[1], xd[2])); 966472b9844SJames Wright else PetscCall(PetscFree3(x[0], x[1], x[2])); 967472b9844SJames Wright 968472b9844SJames Wright { // Reduce CGNS-ordered coordinate nodes to Plex ordering, and set DM's coordinates 969472b9844SJames Wright Vec coord_global; 970472b9844SJames Wright MPI_Datatype unit; 971472b9844SJames Wright PetscScalar *coord_global_array; 972472b9844SJames Wright DM cdm; 973472b9844SJames Wright 974472b9844SJames Wright PetscCall(DMGetCoordinateDM(*dm, &cdm)); 975472b9844SJames Wright PetscCall(DMCreateGlobalVector(cdm, &coord_global)); 976472b9844SJames Wright PetscCall(VecGetArrayWrite(coord_global, &coord_global_array)); 977472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(coordDim, MPIU_SCALAR, &unit)); 978472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 979472b9844SJames Wright PetscCall(PetscSFReduceBegin(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 980472b9844SJames Wright PetscCall(PetscSFReduceEnd(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 981472b9844SJames Wright PetscCall(VecRestoreArrayWrite(coord_global, &coord_global_array)); 982472b9844SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 983472b9844SJames Wright PetscCall(DMSetCoordinates(*dm, coord_global)); 984472b9844SJames Wright PetscCall(VecDestroy(&coord_global)); 985472b9844SJames Wright } 986472b9844SJames Wright PetscCall(PetscFree(coords)); 987472b9844SJames Wright } 988472b9844SJames Wright 989472b9844SJames Wright // -- Set sfNatural for solution vectors in CGNS file 990472b9844SJames Wright // NOTE: We set sfNatural to be the map between the original CGNS ordering of nodes and the Plex ordering of nodes. 991472b9844SJames Wright PetscCall(PetscSFViewFromOptions(plex_to_cgns_sf, NULL, "-sfNatural_init_view")); 992472b9844SJames Wright PetscCall(DMSetNaturalSF(*dm, plex_to_cgns_sf)); 993472b9844SJames Wright PetscCall(DMSetUseNatural(*dm, PETSC_TRUE)); 994472b9844SJames Wright PetscCall(PetscSFDestroy(&plex_to_cgns_sf)); 995472b9844SJames Wright 996472b9844SJames Wright PetscCall(PetscLayoutDestroy(&elem_map)); 997472b9844SJames Wright PetscCall(PetscLayoutDestroy(&vtx_map)); 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9995f34f2dcSJed Brown } 10005f34f2dcSJed Brown 10015f34f2dcSJed Brown // node_l2g must be freed 1002d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g) 1003d71ae5a4SJacob Faibussowitsch { 10045f34f2dcSJed Brown PetscSection local_section; 10055f34f2dcSJed Brown PetscSF point_sf; 10065f34f2dcSJed Brown PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes; 10075f34f2dcSJed Brown PetscMPIInt comm_size; 10085f34f2dcSJed Brown const PetscInt *ilocal, field = 0; 10095f34f2dcSJed Brown 10105f34f2dcSJed Brown PetscFunctionBegin; 10115f34f2dcSJed Brown *num_local_nodes = -1; 10125f34f2dcSJed Brown *num_global_nodes = -1; 10135f34f2dcSJed Brown *nStart = -1; 10145f34f2dcSJed Brown *nEnd = -1; 10155f34f2dcSJed Brown *node_l2g = NULL; 10165f34f2dcSJed Brown 10175f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size)); 10185f34f2dcSJed Brown PetscCall(DMGetLocalSection(dm, &local_section)); 10195f34f2dcSJed Brown PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 10205f34f2dcSJed Brown PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd)); 10215f34f2dcSJed Brown PetscCall(DMGetPointSF(dm, &point_sf)); 10225f34f2dcSJed Brown if (comm_size == 1) nleaves = 0; 10235f34f2dcSJed Brown else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL)); 10245f34f2dcSJed Brown PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp)); 10255f34f2dcSJed Brown 10265f34f2dcSJed Brown PetscInt local_node = 0, owned_node = 0, owned_start = 0; 10275f34f2dcSJed Brown for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 10285f34f2dcSJed Brown PetscInt dof; 10295f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10305f34f2dcSJed Brown PetscAssert(dof % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field dof %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, dof, ncomp); 10315f34f2dcSJed Brown local_node += dof / ncomp; 10325f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 10335f34f2dcSJed Brown leaf++; 10345f34f2dcSJed Brown } else { 10355f34f2dcSJed Brown owned_node += dof / ncomp; 10365f34f2dcSJed Brown } 10375f34f2dcSJed Brown } 10385f34f2dcSJed Brown PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 10395f34f2dcSJed Brown PetscCall(PetscMalloc1(pEnd - pStart, &points)); 10405f34f2dcSJed Brown owned_node = 0; 10415f34f2dcSJed Brown for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 10425f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 10435f34f2dcSJed Brown points[p - pStart] = -1; 10445f34f2dcSJed Brown leaf++; 10455f34f2dcSJed Brown continue; 10465f34f2dcSJed Brown } 10475f34f2dcSJed Brown PetscInt dof, offset; 10485f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10495f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 10505f34f2dcSJed Brown PetscAssert(offset % ncomp == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Field offset %" PetscInt_FMT " must be divisible by components %" PetscInt_FMT, offset, ncomp); 10515f34f2dcSJed Brown points[p - pStart] = owned_start + owned_node; 10525f34f2dcSJed Brown owned_node += dof / ncomp; 10535f34f2dcSJed Brown } 10545f34f2dcSJed Brown if (comm_size > 1) { 10555f34f2dcSJed Brown PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 10565f34f2dcSJed Brown PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 10575f34f2dcSJed Brown } 10585f34f2dcSJed Brown 10595f34f2dcSJed Brown // Set up global indices for each local node 10605f34f2dcSJed Brown PetscCall(PetscMalloc1(local_node, &nodes)); 10615f34f2dcSJed Brown for (PetscInt p = spStart; p < spEnd; p++) { 10625f34f2dcSJed Brown PetscInt dof, offset; 10635f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10645f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 10655f34f2dcSJed Brown for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n; 10665f34f2dcSJed Brown } 10675f34f2dcSJed Brown PetscCall(PetscFree(points)); 10685f34f2dcSJed Brown *num_local_nodes = local_node; 10695f34f2dcSJed Brown *nStart = owned_start; 10705f34f2dcSJed Brown *nEnd = owned_start + owned_node; 1071462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 10725f34f2dcSJed Brown *node_l2g = nodes; 10733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10745f34f2dcSJed Brown } 10755f34f2dcSJed Brown 1076d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer) 1077d71ae5a4SJacob Faibussowitsch { 10785f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 107944698e90SJames Wright PetscInt fvGhostStart; 10805f34f2dcSJed Brown PetscInt topo_dim, coord_dim, num_global_elems; 10815f34f2dcSJed Brown PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd; 10825f34f2dcSJed Brown const PetscInt *node_l2g; 10835f34f2dcSJed Brown Vec coord; 1084b85bf5edSJed Brown DM colloc_dm, cdm; 10855f34f2dcSJed Brown PetscMPIInt size; 10865f34f2dcSJed Brown const char *dm_name; 10875f34f2dcSJed Brown int base, zone; 10885f34f2dcSJed Brown cgsize_t isize[3]; 10895f34f2dcSJed Brown 10905f34f2dcSJed Brown PetscFunctionBegin; 109125760affSJed Brown if (!cgv->file_num) { 109225760affSJed Brown PetscInt time_step; 109325760affSJed Brown PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL)); 109425760affSJed Brown PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step)); 109525760affSJed Brown } 10965f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 10975f34f2dcSJed Brown PetscCall(DMGetDimension(dm, &topo_dim)); 10985f34f2dcSJed Brown PetscCall(DMGetCoordinateDim(dm, &coord_dim)); 10995f34f2dcSJed Brown PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name)); 11004c155f28SJames Wright PetscCallCGNSWrite(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base), dm, viewer); 11015f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, NULL)); 11024c155f28SJames Wright PetscCallCGNSWrite(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)), dm, viewer); 11035f34f2dcSJed Brown 1104b85bf5edSJed Brown { 1105b85bf5edSJed Brown PetscFE fe, fe_coord; 11069bb8f83bSJed Brown PetscClassId ds_id; 1107b85bf5edSJed Brown PetscDualSpace dual_space, dual_space_coord; 11082d1fc720SJed Brown PetscInt num_fields, field_order = -1, field_order_coord; 1109b85bf5edSJed Brown PetscBool is_simplex; 11105b8b2c14SJed Brown PetscCall(DMGetNumFields(dm, &num_fields)); 11119bb8f83bSJed Brown if (num_fields > 0) { 11129bb8f83bSJed Brown PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 11139bb8f83bSJed Brown PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id)); 11142d1fc720SJed Brown if (ds_id != PETSCFE_CLASSID) { 11152d1fc720SJed Brown fe = NULL; 11162d1fc720SJed Brown if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter 11172d1fc720SJed Brown else field_order = 1; // assume vertex-based linear elements 11182d1fc720SJed Brown } 11199bb8f83bSJed Brown } else fe = NULL; 1120b85bf5edSJed Brown if (fe) { 1121b85bf5edSJed Brown PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 1122b85bf5edSJed Brown PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order)); 11232d1fc720SJed Brown } 11245f34f2dcSJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 1125b85bf5edSJed Brown PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord)); 11269bb8f83bSJed Brown { 11279bb8f83bSJed Brown PetscClassId id; 11289bb8f83bSJed Brown PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id)); 11299bb8f83bSJed Brown if (id != PETSCFE_CLASSID) fe_coord = NULL; 11309bb8f83bSJed Brown } 1131b85bf5edSJed Brown if (fe_coord) { 1132b85bf5edSJed Brown PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord)); 1133b85bf5edSJed Brown PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord)); 1134b85bf5edSJed Brown } else field_order_coord = 1; 11352d1fc720SJed Brown if (field_order > 0 && field_order != field_order_coord) { 1136b85bf5edSJed Brown PetscInt quadrature_order = field_order; 1137b85bf5edSJed Brown PetscCall(DMClone(dm, &colloc_dm)); 11387d4eb7abSJed Brown { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied 11391fca310dSJames Wright const PetscSF *face_sfs; 11401fca310dSJames Wright PetscInt num_face_sfs; 11411fca310dSJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs)); 11421fca310dSJames Wright PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs)); 11431fca310dSJames Wright if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal; 11447d4eb7abSJed Brown } 1145b85bf5edSJed Brown PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 1146b85bf5edSJed Brown PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe)); 11474c712d99Sksagiyam PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_FALSE, PETSC_TRUE)); 1148b85bf5edSJed Brown PetscCall(PetscFEDestroy(&fe)); 1149b85bf5edSJed Brown } else { 1150b85bf5edSJed Brown PetscCall(PetscObjectReference((PetscObject)dm)); 1151b85bf5edSJed Brown colloc_dm = dm; 1152b85bf5edSJed Brown } 1153b85bf5edSJed Brown } 1154b85bf5edSJed Brown PetscCall(DMGetCoordinateDM(colloc_dm, &cdm)); 11555f34f2dcSJed Brown PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g)); 1156b85bf5edSJed Brown PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord)); 11575f34f2dcSJed Brown PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 115844698e90SJames Wright PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 115944698e90SJames Wright if (fvGhostStart >= 0) cEnd = fvGhostStart; 11605f34f2dcSJed Brown num_global_elems = cEnd - cStart; 1161462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 11625f34f2dcSJed Brown isize[0] = num_global_nodes; 11635f34f2dcSJed Brown isize[1] = num_global_elems; 11645f34f2dcSJed Brown isize[2] = 0; 11654c155f28SJames Wright PetscCallCGNSWrite(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone), dm, viewer); 11665f34f2dcSJed Brown 11675f34f2dcSJed Brown { 11685f34f2dcSJed Brown const PetscScalar *X; 11695f34f2dcSJed Brown PetscScalar *x; 11705f34f2dcSJed Brown int coord_ids[3]; 11715f34f2dcSJed Brown 11725f34f2dcSJed Brown PetscCall(VecGetArrayRead(coord, &X)); 11735f34f2dcSJed Brown for (PetscInt d = 0; d < coord_dim; d++) { 11745f34f2dcSJed Brown const double exponents[] = {0, 1, 0, 0, 0}; 11755f34f2dcSJed Brown char coord_name[64]; 11765f34f2dcSJed Brown PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d)); 11774c155f28SJames Wright PetscCallCGNSWrite(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]), dm, viewer); 11785f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL)); 11794c155f28SJames Wright PetscCallCGNSWrite(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents), dm, viewer); 11805f34f2dcSJed Brown } 11815f34f2dcSJed Brown 11825f34f2dcSJed Brown int section; 11835f34f2dcSJed Brown cgsize_t e_owned, e_global, e_start, *conn = NULL; 11845f34f2dcSJed Brown const int *perm; 1185e853fb4cSJames Wright CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull); 11865f34f2dcSJed Brown { 11875f34f2dcSJed Brown PetscCall(PetscMalloc1(nEnd - nStart, &x)); 11885f34f2dcSJed Brown for (PetscInt d = 0; d < coord_dim; d++) { 11895f34f2dcSJed Brown for (PetscInt n = 0; n < num_local_nodes; n++) { 11905f34f2dcSJed Brown PetscInt gn = node_l2g[n]; 11915f34f2dcSJed Brown if (gn < nStart || nEnd <= gn) continue; 11925f34f2dcSJed Brown x[gn - nStart] = X[n * coord_dim + d]; 11935f34f2dcSJed Brown } 11945f34f2dcSJed Brown // CGNS nodes use 1-based indexing 11955f34f2dcSJed Brown cgsize_t start = nStart + 1, end = nEnd; 11964c155f28SJames Wright PetscCallCGNSWriteData(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x), dm, viewer); 11975f34f2dcSJed Brown } 11985f34f2dcSJed Brown PetscCall(PetscFree(x)); 11995f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(coord, &X)); 12005f34f2dcSJed Brown } 12015f34f2dcSJed Brown 120259191b1eSJames Wright e_owned = cEnd - cStart; 120344698e90SJames Wright if (e_owned > 0) { 120444698e90SJames Wright DMPolytopeType cell_type; 120544698e90SJames Wright 120644698e90SJames Wright PetscCall(DMPlexGetCellType(dm, cStart, &cell_type)); 12075f34f2dcSJed Brown for (PetscInt i = cStart, c = 0; i < cEnd; i++) { 12085f34f2dcSJed Brown PetscInt closure_dof, *closure_indices, elem_size; 120944698e90SJames Wright 12105f34f2dcSJed Brown PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 12115f34f2dcSJed Brown elem_size = closure_dof / coord_dim; 121244698e90SJames Wright if (!conn) PetscCall(PetscMalloc1(e_owned * elem_size, &conn)); 12135f34f2dcSJed Brown PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm)); 1214ad540459SPierre Jolivet for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1; 12155f34f2dcSJed Brown PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 12165f34f2dcSJed Brown } 121744698e90SJames Wright } 121844698e90SJames Wright 121959191b1eSJames Wright { // Get global element_type (for ranks that do not have owned elements) 122059191b1eSJames Wright PetscInt local_element_type, global_element_type; 122159191b1eSJames Wright 122259191b1eSJames Wright local_element_type = e_owned > 0 ? element_type : -1; 1223462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&local_element_type, &global_element_type, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 122459191b1eSJames Wright 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"); 1225ca328833SJames Wright element_type = (CGNS_ENUMT(ElementType_t))global_element_type; 122659191b1eSJames Wright } 1227462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 12281fb9530eSJeongu Kim 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); 12295f34f2dcSJed Brown e_start = 0; 12301fb9530eSJeongu Kim PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 12314c155f28SJames Wright PetscCallCGNSWrite(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, §ion), dm, viewer); 12324c155f28SJames Wright PetscCallCGNSWriteData(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn), dm, viewer); 12335f34f2dcSJed Brown PetscCall(PetscFree(conn)); 12345f34f2dcSJed Brown 12355f34f2dcSJed Brown cgv->base = base; 12365f34f2dcSJed Brown cgv->zone = zone; 12375f34f2dcSJed Brown cgv->node_l2g = node_l2g; 12385f34f2dcSJed Brown cgv->num_local_nodes = num_local_nodes; 12395f34f2dcSJed Brown cgv->nStart = nStart; 12405f34f2dcSJed Brown cgv->nEnd = nEnd; 12419bb8f83bSJed Brown cgv->eStart = e_start; 12429bb8f83bSJed Brown cgv->eEnd = e_start + e_owned; 12435f34f2dcSJed Brown if (1) { 12445f34f2dcSJed Brown PetscMPIInt rank; 12455f34f2dcSJed Brown int *efield; 12465f34f2dcSJed Brown int sol, field; 12475f34f2dcSJed Brown DMLabel label; 12485f34f2dcSJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 12495f34f2dcSJed Brown PetscCall(PetscMalloc1(e_owned, &efield)); 12505f34f2dcSJed Brown for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank; 12514c155f28SJames Wright PetscCallCGNSWrite(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol), dm, viewer); 12524c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field), dm, viewer); 12535f34f2dcSJed Brown cgsize_t start = e_start + 1, end = e_start + e_owned; 12544c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 12555f34f2dcSJed Brown PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 12565f34f2dcSJed Brown if (label) { 12575f34f2dcSJed Brown for (PetscInt c = cStart; c < cEnd; c++) { 12585f34f2dcSJed Brown PetscInt value; 12595f34f2dcSJed Brown PetscCall(DMLabelGetValue(label, c, &value)); 12605f34f2dcSJed Brown efield[c - cStart] = value; 12615f34f2dcSJed Brown } 12624c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field), dm, viewer); 12634c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 12645f34f2dcSJed Brown } 12655f34f2dcSJed Brown PetscCall(PetscFree(efield)); 12665f34f2dcSJed Brown } 12675f34f2dcSJed Brown } 1268b85bf5edSJed Brown PetscCall(DMDestroy(&colloc_dm)); 12693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12705f34f2dcSJed Brown } 12715f34f2dcSJed Brown 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer) 1273d71ae5a4SJacob Faibussowitsch { 12745f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 12755f34f2dcSJed Brown DM dm; 12765f34f2dcSJed Brown PetscSection section; 127744698e90SJames Wright PetscInt time_step, num_fields, pStart, pEnd, fvGhostStart; 12785f34f2dcSJed Brown PetscReal time, *time_slot; 12799812b6beSJed Brown size_t *step_slot; 12805f34f2dcSJed Brown const PetscScalar *v; 12815f34f2dcSJed Brown char solution_name[PETSC_MAX_PATH_LEN]; 12825f34f2dcSJed Brown int sol; 12835f34f2dcSJed Brown 12845f34f2dcSJed Brown PetscFunctionBegin; 12855f34f2dcSJed Brown PetscCall(VecGetDM(V, &dm)); 1286ec2db9e4SJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 1287ec2db9e4SJames Wright PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 128844698e90SJames Wright PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 128944698e90SJames Wright if (fvGhostStart >= 0) pEnd = fvGhostStart; 1290ec2db9e4SJames Wright 12915f34f2dcSJed Brown if (!cgv->node_l2g) PetscCall(DMView(dm, viewer)); 1292ec2db9e4SJames Wright if (!cgv->grid_loc) { // Determine if writing to cell-centers or to nodes 1293ec2db9e4SJames Wright PetscInt cStart, cEnd; 1294ec2db9e4SJames Wright PetscInt local_grid_loc, global_grid_loc; 1295ec2db9e4SJames Wright 1296ec2db9e4SJames Wright PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 129744698e90SJames Wright if (fvGhostStart >= 0) cEnd = fvGhostStart; 1298ec2db9e4SJames Wright if (cgv->num_local_nodes == 0) local_grid_loc = -1; 1299ec2db9e4SJames Wright else if (cStart == pStart && cEnd == pEnd) local_grid_loc = CGNS_ENUMV(CellCenter); 1300ec2db9e4SJames Wright else local_grid_loc = CGNS_ENUMV(Vertex); 1301ec2db9e4SJames Wright 1302462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&local_grid_loc, &global_grid_loc, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 1303ec2db9e4SJames Wright if (local_grid_loc != -1) 1304ec2db9e4SJames Wright 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); 1305ec2db9e4SJames Wright 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); 1306ca328833SJames Wright cgv->grid_loc = (CGNS_ENUMT(GridLocation_t))global_grid_loc; 1307ec2db9e4SJames Wright } 1308ec2db9e4SJames Wright if (!cgv->nodal_field) { 1309ec2db9e4SJames Wright switch (cgv->grid_loc) { 1310ec2db9e4SJames Wright case CGNS_ENUMV(Vertex): { 1311ec2db9e4SJames Wright PetscCall(PetscMalloc1(cgv->nEnd - cgv->nStart, &cgv->nodal_field)); 1312ec2db9e4SJames Wright } break; 1313ec2db9e4SJames Wright case CGNS_ENUMV(CellCenter): { 1314ec2db9e4SJames Wright PetscCall(PetscMalloc1(cgv->eEnd - cgv->eStart, &cgv->nodal_field)); 1315ec2db9e4SJames Wright } break; 1316ec2db9e4SJames Wright default: 1317ec2db9e4SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 1318ec2db9e4SJames Wright } 1319ec2db9e4SJames Wright } 13205f34f2dcSJed Brown if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times)); 13219812b6beSJed Brown if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps)); 13225f34f2dcSJed Brown 13235f34f2dcSJed Brown PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time)); 13249371c9d4SSatish Balay if (time_step < 0) { 13259371c9d4SSatish Balay time_step = 0; 13269371c9d4SSatish Balay time = 0.; 13279371c9d4SSatish Balay } 1328*dfa0de3dSJames Wright // Avoid "Duplicate child name found" error by not writing an already-written solution. 1329*dfa0de3dSJames Wright // This usually occurs when a solution is written and then diverges on the very next timestep. 1330*dfa0de3dSJames Wright if (time_step == cgv->previous_output_step) PetscFunctionReturn(PETSC_SUCCESS); 1331*dfa0de3dSJames Wright 13325f34f2dcSJed Brown PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot)); 13335f34f2dcSJed Brown *time_slot = time; 13349812b6beSJed Brown PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot)); 1335*dfa0de3dSJames Wright *step_slot = cgv->previous_output_step = time_step; 13365f34f2dcSJed Brown PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step)); 13374c155f28SJames Wright PetscCallCGNSWrite(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, cgv->grid_loc, &sol), V, viewer); 13385f34f2dcSJed Brown PetscCall(VecGetArrayRead(V, &v)); 133900a86070SJed Brown PetscCall(PetscSectionGetNumFields(section, &num_fields)); 134000a86070SJed Brown for (PetscInt field = 0; field < num_fields; field++) { 134100a86070SJed Brown PetscInt ncomp; 13429bb8f83bSJed Brown const char *field_name; 13439bb8f83bSJed Brown PetscCall(PetscSectionGetFieldName(section, field, &field_name)); 134400a86070SJed Brown PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp)); 13455f34f2dcSJed Brown for (PetscInt comp = 0; comp < ncomp; comp++) { 13465f34f2dcSJed Brown int cgfield; 13475f34f2dcSJed Brown const char *comp_name; 13489bb8f83bSJed Brown char cgns_field_name[32]; // CGNS max field name is 32 13495f34f2dcSJed Brown CGNS_ENUMT(DataType_t) datatype; 13505f34f2dcSJed Brown PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name)); 135144698e90SJames Wright 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)); 13529bb8f83bSJed Brown else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name)); 13539bb8f83bSJed Brown else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name)); 13545f34f2dcSJed Brown PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype)); 13554c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield), V, viewer); 135600a86070SJed Brown for (PetscInt p = pStart, n = 0; p < pEnd; p++) { 135700a86070SJed Brown PetscInt off, dof; 135800a86070SJed Brown PetscCall(PetscSectionGetFieldDof(section, p, field, &dof)); 135900a86070SJed Brown if (dof == 0) continue; 136000a86070SJed Brown PetscCall(PetscSectionGetFieldOffset(section, p, field, &off)); 136100a86070SJed Brown for (PetscInt c = comp; c < dof; c += ncomp, n++) { 1362ec2db9e4SJames Wright switch (cgv->grid_loc) { 13639bb8f83bSJed Brown case CGNS_ENUMV(Vertex): { 13645f34f2dcSJed Brown PetscInt gn = cgv->node_l2g[n]; 13655f34f2dcSJed Brown if (gn < cgv->nStart || cgv->nEnd <= gn) continue; 136600a86070SJed Brown cgv->nodal_field[gn - cgv->nStart] = v[off + c]; 13679bb8f83bSJed Brown } break; 13689bb8f83bSJed Brown case CGNS_ENUMV(CellCenter): { 13699bb8f83bSJed Brown cgv->nodal_field[n] = v[off + c]; 13709bb8f83bSJed Brown } break; 13719bb8f83bSJed Brown default: 13729bb8f83bSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations"); 13739bb8f83bSJed Brown } 137400a86070SJed Brown } 13755f34f2dcSJed Brown } 13765f34f2dcSJed Brown // CGNS nodes use 1-based indexing 1377ec2db9e4SJames Wright cgsize_t start, end; 1378ec2db9e4SJames Wright switch (cgv->grid_loc) { 1379ec2db9e4SJames Wright case CGNS_ENUMV(Vertex): { 1380ec2db9e4SJames Wright start = cgv->nStart + 1; 1381ec2db9e4SJames Wright end = cgv->nEnd; 1382ec2db9e4SJames Wright } break; 1383ec2db9e4SJames Wright case CGNS_ENUMV(CellCenter): { 13849bb8f83bSJed Brown start = cgv->eStart + 1; 13859bb8f83bSJed Brown end = cgv->eEnd; 1386ec2db9e4SJames Wright } break; 1387ec2db9e4SJames Wright default: 1388ec2db9e4SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 13899bb8f83bSJed Brown } 13904c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field), V, viewer); 13915f34f2dcSJed Brown } 139200a86070SJed Brown } 13935f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(V, &v)); 139425760affSJed Brown PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer)); 13953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13965f34f2dcSJed Brown } 1397472b9844SJames Wright 1398472b9844SJames Wright PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec V, PetscViewer viewer) 1399472b9844SJames Wright { 1400472b9844SJames Wright MPI_Comm comm; 1401472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 1402472b9844SJames Wright PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 1403472b9844SJames Wright int cgid = cgv->file_num; 1404472b9844SJames Wright PetscBool use_parallel_viewer = PETSC_FALSE; 1405472b9844SJames Wright int z = 1; // Only support one zone 1406472b9844SJames Wright int B = 1; // Only support one base 1407472b9844SJames Wright int numComp; 1408472b9844SJames Wright PetscInt V_numComps, mystartv, myendv, myownedv; 1409472b9844SJames Wright 1410472b9844SJames Wright PetscFunctionBeginUser; 1411472b9844SJames Wright PetscCall(PetscObjectGetComm((PetscObject)V, &comm)); 1412472b9844SJames Wright 1413472b9844SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 1414472b9844SJames Wright 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"); 1415472b9844SJames Wright 1416472b9844SJames Wright { // Get CGNS node ownership information 1417472b9844SJames Wright int nbases, nzones; 1418472b9844SJames Wright PetscInt NVertices; 1419472b9844SJames Wright PetscLayout vtx_map; 1420472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 1421472b9844SJames Wright 14224c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), V, viewer); 1423472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 14244c155f28SJames Wright PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), V, viewer); 1425472b9844SJames Wright PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "limited to one zone %d", (int)nzones); 1426472b9844SJames Wright 14274c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), V, viewer); 1428472b9844SJames Wright NVertices = sizes[0]; 1429472b9844SJames Wright 1430472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 1431472b9844SJames Wright PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 1432472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 1433472b9844SJames Wright PetscCall(PetscLayoutDestroy(&vtx_map)); 1434472b9844SJames Wright } 1435472b9844SJames Wright 1436472b9844SJames Wright { // -- Read data from file into Vec 1437472b9844SJames Wright PetscScalar *fields = NULL; 1438472b9844SJames Wright PetscSF sfNatural; 1439472b9844SJames Wright 1440472b9844SJames Wright { // Check compatibility between sfNatural and the data source and sink 1441472b9844SJames Wright DM dm; 1442472b9844SJames Wright PetscInt nleaves, nroots, V_local_size; 1443472b9844SJames Wright 1444472b9844SJames Wright PetscCall(VecGetDM(V, &dm)); 1445472b9844SJames Wright PetscCall(DMGetNaturalSF(dm, &sfNatural)); 1446472b9844SJames Wright PetscCheck(sfNatural, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM of Vec must have sfNatural"); 1447472b9844SJames Wright PetscCall(PetscSFGetGraph(sfNatural, &nroots, &nleaves, NULL, NULL)); 1448472b9844SJames Wright PetscCall(VecGetLocalSize(V, &V_local_size)); 1449472b9844SJames Wright 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); 1450472b9844SJames Wright 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); 1451472b9844SJames Wright V_numComps = V_local_size / nroots; 1452472b9844SJames Wright } 1453472b9844SJames Wright 1454472b9844SJames Wright { // Read data into component-major ordering 1455472b9844SJames Wright int isol, numSols; 1456472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 1457472b9844SJames Wright double *fields_CGNS; 1458472b9844SJames Wright 14594c155f28SJames Wright PetscCallCGNSRead(cg_nsols(cgid, B, z, &numSols), V, viewer); 1460472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &isol)); 14614c155f28SJames Wright PetscCallCGNSRead(cg_nfields(cgid, B, z, isol, &numComp), V, viewer); 1462472b9844SJames Wright 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); 1463472b9844SJames Wright 1464472b9844SJames Wright cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 1465472b9844SJames Wright cgsize_t range_max[3] = {myendv, 1, 1}; 1466472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * numComp, &fields_CGNS)); 1467472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * numComp, &fields)); 1468472b9844SJames Wright for (int d = 0; d < numComp; ++d) { 14694c155f28SJames Wright PetscCallCGNSRead(cg_field_info(cgid, B, z, isol, (d + 1), &datatype, buffer), V, viewer); 1470472b9844SJames Wright PetscCheck(datatype == CGNS_ENUMV(RealDouble), PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Field %s in file is not of type double", buffer); 14714c155f28SJames Wright PetscCallCGNSReadData(cgp_field_read_data(cgid, B, z, isol, (d + 1), range_min, range_max, &fields_CGNS[d * myownedv]), V, viewer); 1472472b9844SJames Wright } 1473472b9844SJames Wright for (int d = 0; d < numComp; ++d) { 1474472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) fields[v * numComp + d] = fields_CGNS[d * myownedv + v]; 1475472b9844SJames Wright } 1476472b9844SJames Wright PetscCall(PetscFree(fields_CGNS)); 1477472b9844SJames Wright } 1478472b9844SJames Wright 1479472b9844SJames Wright { // Reduce fields into Vec array 1480472b9844SJames Wright PetscScalar *V_array; 1481472b9844SJames Wright MPI_Datatype fieldtype; 1482472b9844SJames Wright 1483472b9844SJames Wright PetscCall(VecGetArrayWrite(V, &V_array)); 1484472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(numComp, MPIU_SCALAR, &fieldtype)); 1485472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&fieldtype)); 1486472b9844SJames Wright PetscCall(PetscSFReduceBegin(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1487472b9844SJames Wright PetscCall(PetscSFReduceEnd(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1488472b9844SJames Wright PetscCallMPI(MPI_Type_free(&fieldtype)); 1489472b9844SJames Wright PetscCall(VecRestoreArrayWrite(V, &V_array)); 1490472b9844SJames Wright } 1491472b9844SJames Wright PetscCall(PetscFree(fields)); 1492472b9844SJames Wright } 1493472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1494472b9844SJames Wright } 1495