15f34f2dcSJed Brown #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 25f34f2dcSJed Brown #include <petsc/private/viewercgnsimpl.h> 35f34f2dcSJed Brown 45f34f2dcSJed Brown #include <pcgnslib.h> 55f34f2dcSJed Brown #include <cgns_io.h> 65f34f2dcSJed Brown 75f34f2dcSJed Brown #if !defined(CGNS_ENUMT) 85f34f2dcSJed Brown #define CGNS_ENUMT(a) a 95f34f2dcSJed Brown #endif 105f34f2dcSJed Brown #if !defined(CGNS_ENUMV) 115f34f2dcSJed Brown #define CGNS_ENUMV(a) a 125f34f2dcSJed Brown #endif 135f34f2dcSJed Brown // Permute plex closure ordering to CGNS 14d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCGNSGetPermutation_Internal(DMPolytopeType cell_type, PetscInt closure_size, CGNS_ENUMT(ElementType_t) * element_type, const int **perm) 15d71ae5a4SJacob Faibussowitsch { 16472b9844SJames Wright CGNS_ENUMT(ElementType_t) element_type_tmp; 17472b9844SJames Wright 18472b9844SJames Wright // https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 195f34f2dcSJed Brown static const int bar_2[2] = {0, 1}; 205f34f2dcSJed Brown static const int bar_3[3] = {1, 2, 0}; 215f34f2dcSJed Brown static const int bar_4[4] = {2, 3, 0, 1}; 225f34f2dcSJed Brown static const int bar_5[5] = {3, 4, 0, 1, 2}; 235f34f2dcSJed Brown static const int tri_3[3] = {0, 1, 2}; 245f34f2dcSJed Brown static const int tri_6[6] = {3, 4, 5, 0, 1, 2}; 255f34f2dcSJed Brown static const int tri_10[10] = {7, 8, 9, 1, 2, 3, 4, 5, 6, 0}; 265f34f2dcSJed Brown static const int quad_4[4] = {0, 1, 2, 3}; 275f34f2dcSJed Brown static const int quad_9[9] = { 285f34f2dcSJed Brown 5, 6, 7, 8, // vertices 295f34f2dcSJed Brown 1, 2, 3, 4, // edges 305f34f2dcSJed Brown 0, // center 315f34f2dcSJed Brown }; 325f34f2dcSJed Brown static const int quad_16[] = { 335f34f2dcSJed Brown 12, 13, 14, 15, // vertices 345f34f2dcSJed Brown 4, 5, 6, 7, 8, 9, 10, 11, // edges 355f34f2dcSJed Brown 0, 1, 3, 2, // centers 365f34f2dcSJed Brown }; 375f34f2dcSJed Brown static const int quad_25[] = { 385f34f2dcSJed Brown 21, 22, 23, 24, // vertices 395f34f2dcSJed Brown 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, // edges 405f34f2dcSJed Brown 0, 1, 2, 5, 8, 7, 6, 3, 4, // centers 415f34f2dcSJed Brown }; 425f34f2dcSJed Brown static const int tetra_4[4] = {0, 2, 1, 3}; 435f34f2dcSJed Brown static const int tetra_10[10] = {6, 8, 7, 9, 2, 1, 0, 3, 5, 4}; 445f34f2dcSJed Brown static const int tetra_20[20] = { 455f34f2dcSJed Brown 16, 18, 17, 19, // vertices 465f34f2dcSJed Brown 9, 8, 7, 6, 5, 4, // bottom edges 475f34f2dcSJed Brown 10, 11, 14, 15, 13, 12, // side edges 485f34f2dcSJed Brown 0, 2, 3, 1, // faces 495f34f2dcSJed Brown }; 505f34f2dcSJed Brown static const int hexa_8[8] = {0, 3, 2, 1, 4, 5, 6, 7}; 515f34f2dcSJed Brown static const int hexa_27[27] = { 525f34f2dcSJed Brown 19, 22, 21, 20, 23, 24, 25, 26, // vertices 535f34f2dcSJed Brown 10, 9, 8, 7, // bottom edges 545f34f2dcSJed Brown 16, 15, 18, 17, // mid edges 555f34f2dcSJed Brown 11, 12, 13, 14, // top edges 565f34f2dcSJed Brown 1, 3, 5, 4, 6, 2, // faces 575f34f2dcSJed Brown 0, // center 585f34f2dcSJed Brown }; 599371c9d4SSatish Balay static const int hexa_64[64] = { 609371c9d4SSatish 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 615f34f2dcSJed Brown 56, 59, 58, 57, 60, 61, 62, 63, // vertices 625f34f2dcSJed Brown 39, 38, 37, 36, 35, 34, 33, 32, // bottom edges 635f34f2dcSJed Brown 51, 50, 48, 49, 52, 53, 55, 54, // mid edges; Paraview needs edge 21-22 swapped with 23-24 645f34f2dcSJed Brown 40, 41, 42, 43, 44, 45, 46, 47, // top edges 655f34f2dcSJed Brown 8, 10, 11, 9, // z-minus face 665f34f2dcSJed Brown 16, 17, 19, 18, // y-minus face 675f34f2dcSJed Brown 24, 25, 27, 26, // x-plus face 685f34f2dcSJed Brown 20, 21, 23, 22, // y-plus face 695f34f2dcSJed Brown 30, 28, 29, 31, // x-minus face 705f34f2dcSJed Brown 12, 13, 15, 14, // z-plus face 715f34f2dcSJed Brown 0, 1, 3, 2, 4, 5, 7, 6, // center 725f34f2dcSJed Brown }; 735f34f2dcSJed Brown 745f34f2dcSJed Brown PetscFunctionBegin; 75472b9844SJames Wright element_type_tmp = CGNS_ENUMV(ElementTypeNull); 765f34f2dcSJed Brown *perm = NULL; 775f34f2dcSJed Brown switch (cell_type) { 785f34f2dcSJed Brown case DM_POLYTOPE_SEGMENT: 795f34f2dcSJed Brown switch (closure_size) { 80d71ae5a4SJacob Faibussowitsch case 2: 81472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_2); 82d71ae5a4SJacob Faibussowitsch *perm = bar_2; 83b30057acSJames Wright break; 84d71ae5a4SJacob Faibussowitsch case 3: 85472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_3); 86d71ae5a4SJacob Faibussowitsch *perm = bar_3; 87b30057acSJames Wright break; 885f34f2dcSJed Brown case 4: 89472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_4); 905f34f2dcSJed Brown *perm = bar_4; 915f34f2dcSJed Brown break; 925f34f2dcSJed Brown case 5: 93472b9844SJames Wright element_type_tmp = CGNS_ENUMV(BAR_5); 945f34f2dcSJed Brown *perm = bar_5; 955f34f2dcSJed Brown break; 96d71ae5a4SJacob Faibussowitsch default: 97d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 985f34f2dcSJed Brown } 995f34f2dcSJed Brown break; 1005f34f2dcSJed Brown case DM_POLYTOPE_TRIANGLE: 1015f34f2dcSJed Brown switch (closure_size) { 1025f34f2dcSJed Brown case 3: 103472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_3); 1045f34f2dcSJed Brown *perm = tri_3; 1055f34f2dcSJed Brown break; 1065f34f2dcSJed Brown case 6: 107472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_6); 1085f34f2dcSJed Brown *perm = tri_6; 1095f34f2dcSJed Brown break; 1105f34f2dcSJed Brown case 10: 111472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TRI_10); 1125f34f2dcSJed Brown *perm = tri_10; 1135f34f2dcSJed Brown break; 114d71ae5a4SJacob Faibussowitsch default: 115d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1165f34f2dcSJed Brown } 1175f34f2dcSJed Brown break; 1185f34f2dcSJed Brown case DM_POLYTOPE_QUADRILATERAL: 1195f34f2dcSJed Brown switch (closure_size) { 1205f34f2dcSJed Brown case 4: 121472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_4); 1225f34f2dcSJed Brown *perm = quad_4; 1235f34f2dcSJed Brown break; 1245f34f2dcSJed Brown case 9: 125472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_9); 1265f34f2dcSJed Brown *perm = quad_9; 1275f34f2dcSJed Brown break; 1285f34f2dcSJed Brown case 16: 129472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_16); 1305f34f2dcSJed Brown *perm = quad_16; 1315f34f2dcSJed Brown break; 1325f34f2dcSJed Brown case 25: 133472b9844SJames Wright element_type_tmp = CGNS_ENUMV(QUAD_25); 1345f34f2dcSJed Brown *perm = quad_25; 1355f34f2dcSJed Brown break; 136d71ae5a4SJacob Faibussowitsch default: 137d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1385f34f2dcSJed Brown } 1395f34f2dcSJed Brown break; 1405f34f2dcSJed Brown case DM_POLYTOPE_TETRAHEDRON: 1415f34f2dcSJed Brown switch (closure_size) { 1425f34f2dcSJed Brown case 4: 143472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_4); 1445f34f2dcSJed Brown *perm = tetra_4; 1455f34f2dcSJed Brown break; 1465f34f2dcSJed Brown case 10: 147472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_10); 1485f34f2dcSJed Brown *perm = tetra_10; 1495f34f2dcSJed Brown break; 1505f34f2dcSJed Brown case 20: 151472b9844SJames Wright element_type_tmp = CGNS_ENUMV(TETRA_20); 1525f34f2dcSJed Brown *perm = tetra_20; 1535f34f2dcSJed Brown break; 154d71ae5a4SJacob Faibussowitsch default: 155d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1565f34f2dcSJed Brown } 1575f34f2dcSJed Brown break; 1585f34f2dcSJed Brown case DM_POLYTOPE_HEXAHEDRON: 1595f34f2dcSJed Brown switch (closure_size) { 1605f34f2dcSJed Brown case 8: 161472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_8); 1625f34f2dcSJed Brown *perm = hexa_8; 1635f34f2dcSJed Brown break; 1645f34f2dcSJed Brown case 27: 165472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_27); 1665f34f2dcSJed Brown *perm = hexa_27; 1675f34f2dcSJed Brown break; 1685f34f2dcSJed Brown case 64: 169472b9844SJames Wright element_type_tmp = CGNS_ENUMV(HEXA_64); 1705f34f2dcSJed Brown *perm = hexa_64; 1715f34f2dcSJed Brown break; 172d71ae5a4SJacob Faibussowitsch default: 173d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1745f34f2dcSJed Brown } 1755f34f2dcSJed Brown break; 176d71ae5a4SJacob Faibussowitsch default: 177d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cell type %s with closure size %" PetscInt_FMT, DMPolytopeTypes[cell_type], closure_size); 1785f34f2dcSJed Brown } 179472b9844SJames Wright if (element_type) *element_type = element_type_tmp; 180472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 181472b9844SJames Wright } 182472b9844SJames Wright 183472b9844SJames Wright /* 184472b9844SJames Wright Input Parameters: 185472b9844SJames Wright + cellType - The CGNS-defined element type 186472b9844SJames Wright 187472b9844SJames Wright Output Parameters: 188472b9844SJames Wright + dmcelltype - The equivalent DMPolytopeType for the cellType 189472b9844SJames Wright . numCorners - Number of corners of the polytope 190472b9844SJames Wright - dim - The topological dimension of the polytope 191472b9844SJames Wright 192472b9844SJames Wright CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 193472b9844SJames Wright */ 194472b9844SJames Wright static inline PetscErrorCode CGNSElementTypeGetTopologyInfo(CGNS_ENUMT(ElementType_t) cellType, DMPolytopeType *dmcelltype, PetscInt *numCorners, PetscInt *dim) 195472b9844SJames Wright { 196472b9844SJames Wright DMPolytopeType _dmcelltype; 197472b9844SJames Wright 198b30057acSJames Wright PetscFunctionBegin; 199472b9844SJames Wright switch (cellType) { 200472b9844SJames Wright case CGNS_ENUMV(BAR_2): 201472b9844SJames Wright case CGNS_ENUMV(BAR_3): 202472b9844SJames Wright case CGNS_ENUMV(BAR_4): 203472b9844SJames Wright case CGNS_ENUMV(BAR_5): 204472b9844SJames Wright _dmcelltype = DM_POLYTOPE_SEGMENT; 205472b9844SJames Wright break; 206472b9844SJames Wright case CGNS_ENUMV(TRI_3): 207472b9844SJames Wright case CGNS_ENUMV(TRI_6): 208472b9844SJames Wright case CGNS_ENUMV(TRI_9): 209472b9844SJames Wright case CGNS_ENUMV(TRI_10): 210472b9844SJames Wright case CGNS_ENUMV(TRI_12): 211472b9844SJames Wright case CGNS_ENUMV(TRI_15): 212472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TRIANGLE; 213472b9844SJames Wright break; 214472b9844SJames Wright case CGNS_ENUMV(QUAD_4): 215472b9844SJames Wright case CGNS_ENUMV(QUAD_8): 216472b9844SJames Wright case CGNS_ENUMV(QUAD_9): 217472b9844SJames Wright case CGNS_ENUMV(QUAD_12): 218472b9844SJames Wright case CGNS_ENUMV(QUAD_16): 219472b9844SJames Wright case CGNS_ENUMV(QUAD_P4_16): 220472b9844SJames Wright case CGNS_ENUMV(QUAD_25): 221472b9844SJames Wright _dmcelltype = DM_POLYTOPE_QUADRILATERAL; 222472b9844SJames Wright break; 223472b9844SJames Wright case CGNS_ENUMV(TETRA_4): 224472b9844SJames Wright case CGNS_ENUMV(TETRA_10): 225472b9844SJames Wright case CGNS_ENUMV(TETRA_16): 226472b9844SJames Wright case CGNS_ENUMV(TETRA_20): 227472b9844SJames Wright case CGNS_ENUMV(TETRA_22): 228472b9844SJames Wright case CGNS_ENUMV(TETRA_34): 229472b9844SJames Wright case CGNS_ENUMV(TETRA_35): 230472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TETRAHEDRON; 231472b9844SJames Wright break; 232472b9844SJames Wright case CGNS_ENUMV(PYRA_5): 233472b9844SJames Wright case CGNS_ENUMV(PYRA_13): 234472b9844SJames Wright case CGNS_ENUMV(PYRA_14): 235472b9844SJames Wright case CGNS_ENUMV(PYRA_21): 236472b9844SJames Wright case CGNS_ENUMV(PYRA_29): 237472b9844SJames Wright case CGNS_ENUMV(PYRA_P4_29): 238472b9844SJames Wright case CGNS_ENUMV(PYRA_30): 239472b9844SJames Wright case CGNS_ENUMV(PYRA_50): 240472b9844SJames Wright case CGNS_ENUMV(PYRA_55): 241472b9844SJames Wright _dmcelltype = DM_POLYTOPE_PYRAMID; 242472b9844SJames Wright break; 243472b9844SJames Wright case CGNS_ENUMV(PENTA_6): 244472b9844SJames Wright case CGNS_ENUMV(PENTA_15): 245472b9844SJames Wright case CGNS_ENUMV(PENTA_18): 246472b9844SJames Wright case CGNS_ENUMV(PENTA_24): 247472b9844SJames Wright case CGNS_ENUMV(PENTA_33): 248472b9844SJames Wright case CGNS_ENUMV(PENTA_38): 249472b9844SJames Wright case CGNS_ENUMV(PENTA_40): 250472b9844SJames Wright case CGNS_ENUMV(PENTA_66): 251472b9844SJames Wright case CGNS_ENUMV(PENTA_75): 252472b9844SJames Wright _dmcelltype = DM_POLYTOPE_TRI_PRISM; 253472b9844SJames Wright break; 254472b9844SJames Wright case CGNS_ENUMV(HEXA_8): 255472b9844SJames Wright case CGNS_ENUMV(HEXA_20): 256472b9844SJames Wright case CGNS_ENUMV(HEXA_27): 257472b9844SJames Wright case CGNS_ENUMV(HEXA_32): 258472b9844SJames Wright case CGNS_ENUMV(HEXA_44): 259472b9844SJames Wright case CGNS_ENUMV(HEXA_56): 260472b9844SJames Wright case CGNS_ENUMV(HEXA_64): 261472b9844SJames Wright case CGNS_ENUMV(HEXA_98): 262472b9844SJames Wright case CGNS_ENUMV(HEXA_125): 263472b9844SJames Wright _dmcelltype = DM_POLYTOPE_HEXAHEDRON; 264472b9844SJames Wright break; 265472b9844SJames Wright case CGNS_ENUMV(MIXED): 266472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 267472b9844SJames Wright default: 268472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: %d", (int)cellType); 269472b9844SJames Wright } 270472b9844SJames Wright 271472b9844SJames Wright if (dmcelltype) *dmcelltype = _dmcelltype; 272472b9844SJames Wright if (numCorners) *numCorners = DMPolytopeTypeGetNumVertices(_dmcelltype); 273472b9844SJames Wright if (dim) *dim = DMPolytopeTypeGetDim(_dmcelltype); 274472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 275472b9844SJames Wright } 276472b9844SJames Wright 277472b9844SJames Wright /* 278472b9844SJames Wright Input Parameters: 279472b9844SJames Wright + cellType - The CGNS-defined cell type 280472b9844SJames Wright 281472b9844SJames Wright Output Parameters: 282472b9844SJames Wright + numClosure - Number of nodes that define the function space on the cell 283472b9844SJames Wright - pOrder - The polynomial order of the cell 284472b9844SJames Wright 285472b9844SJames Wright CGNS elements defined in: https://cgns.github.io/CGNS_docs_current/sids/conv.html#unstructgrid 286472b9844SJames Wright 287472b9844SJames Wright Note: we only support "full" elements, ie. not seredipity elements 288472b9844SJames Wright */ 289472b9844SJames Wright static inline PetscErrorCode CGNSElementTypeGetDiscretizationInfo(CGNS_ENUMT(ElementType_t) cellType, PetscInt *numClosure, PetscInt *pOrder) 290472b9844SJames Wright { 291472b9844SJames Wright PetscInt _numClosure, _pOrder; 292472b9844SJames Wright 293b30057acSJames Wright PetscFunctionBegin; 294472b9844SJames Wright switch (cellType) { 295472b9844SJames Wright case CGNS_ENUMV(BAR_2): 296472b9844SJames Wright _numClosure = 2; 297472b9844SJames Wright _pOrder = 1; 298472b9844SJames Wright break; 299472b9844SJames Wright case CGNS_ENUMV(BAR_3): 300472b9844SJames Wright _numClosure = 3; 301472b9844SJames Wright _pOrder = 2; 302472b9844SJames Wright break; 303472b9844SJames Wright case CGNS_ENUMV(BAR_4): 304472b9844SJames Wright _numClosure = 4; 305472b9844SJames Wright _pOrder = 3; 306472b9844SJames Wright break; 307472b9844SJames Wright case CGNS_ENUMV(BAR_5): 308472b9844SJames Wright _numClosure = 5; 309472b9844SJames Wright _pOrder = 4; 310472b9844SJames Wright break; 311472b9844SJames Wright case CGNS_ENUMV(TRI_3): 312472b9844SJames Wright _numClosure = 3; 313472b9844SJames Wright _pOrder = 1; 314472b9844SJames Wright break; 315472b9844SJames Wright case CGNS_ENUMV(TRI_6): 316472b9844SJames Wright _numClosure = 6; 317472b9844SJames Wright _pOrder = 2; 318472b9844SJames Wright break; 319472b9844SJames Wright case CGNS_ENUMV(TRI_10): 320472b9844SJames Wright _numClosure = 10; 321472b9844SJames Wright _pOrder = 3; 322472b9844SJames Wright break; 323472b9844SJames Wright case CGNS_ENUMV(TRI_15): 324472b9844SJames Wright _numClosure = 15; 325472b9844SJames Wright _pOrder = 4; 326472b9844SJames Wright break; 327472b9844SJames Wright case CGNS_ENUMV(QUAD_4): 328472b9844SJames Wright _numClosure = 4; 329472b9844SJames Wright _pOrder = 1; 330472b9844SJames Wright break; 331472b9844SJames Wright case CGNS_ENUMV(QUAD_9): 332472b9844SJames Wright _numClosure = 9; 333472b9844SJames Wright _pOrder = 2; 334472b9844SJames Wright break; 335472b9844SJames Wright case CGNS_ENUMV(QUAD_16): 336472b9844SJames Wright _numClosure = 16; 337472b9844SJames Wright _pOrder = 3; 338472b9844SJames Wright break; 339472b9844SJames Wright case CGNS_ENUMV(QUAD_25): 340472b9844SJames Wright _numClosure = 25; 341472b9844SJames Wright _pOrder = 4; 342472b9844SJames Wright break; 343472b9844SJames Wright case CGNS_ENUMV(TETRA_4): 344472b9844SJames Wright _numClosure = 4; 345472b9844SJames Wright _pOrder = 1; 346472b9844SJames Wright break; 347472b9844SJames Wright case CGNS_ENUMV(TETRA_10): 348472b9844SJames Wright _numClosure = 10; 349472b9844SJames Wright _pOrder = 2; 350472b9844SJames Wright break; 351472b9844SJames Wright case CGNS_ENUMV(TETRA_20): 352472b9844SJames Wright _numClosure = 20; 353472b9844SJames Wright _pOrder = 3; 354472b9844SJames Wright break; 355472b9844SJames Wright case CGNS_ENUMV(TETRA_35): 356472b9844SJames Wright _numClosure = 35; 357472b9844SJames Wright _pOrder = 4; 358472b9844SJames Wright break; 359472b9844SJames Wright case CGNS_ENUMV(PYRA_5): 360472b9844SJames Wright _numClosure = 5; 361472b9844SJames Wright _pOrder = 1; 362472b9844SJames Wright break; 363472b9844SJames Wright case CGNS_ENUMV(PYRA_14): 364472b9844SJames Wright _numClosure = 14; 365472b9844SJames Wright _pOrder = 2; 366472b9844SJames Wright break; 367472b9844SJames Wright case CGNS_ENUMV(PYRA_30): 368472b9844SJames Wright _numClosure = 30; 369472b9844SJames Wright _pOrder = 3; 370472b9844SJames Wright break; 371472b9844SJames Wright case CGNS_ENUMV(PYRA_55): 372472b9844SJames Wright _numClosure = 55; 373472b9844SJames Wright _pOrder = 4; 374472b9844SJames Wright break; 375472b9844SJames Wright case CGNS_ENUMV(PENTA_6): 376472b9844SJames Wright _numClosure = 6; 377472b9844SJames Wright _pOrder = 1; 378472b9844SJames Wright break; 379472b9844SJames Wright case CGNS_ENUMV(PENTA_18): 380472b9844SJames Wright _numClosure = 18; 381472b9844SJames Wright _pOrder = 2; 382472b9844SJames Wright break; 383472b9844SJames Wright case CGNS_ENUMV(PENTA_40): 384472b9844SJames Wright _numClosure = 40; 385472b9844SJames Wright _pOrder = 3; 386472b9844SJames Wright break; 387472b9844SJames Wright case CGNS_ENUMV(PENTA_75): 388472b9844SJames Wright _numClosure = 75; 389472b9844SJames Wright _pOrder = 4; 390472b9844SJames Wright break; 391472b9844SJames Wright case CGNS_ENUMV(HEXA_8): 392472b9844SJames Wright _numClosure = 8; 393472b9844SJames Wright _pOrder = 1; 394472b9844SJames Wright break; 395472b9844SJames Wright case CGNS_ENUMV(HEXA_27): 396472b9844SJames Wright _numClosure = 27; 397472b9844SJames Wright _pOrder = 2; 398472b9844SJames Wright break; 399472b9844SJames Wright case CGNS_ENUMV(HEXA_64): 400472b9844SJames Wright _numClosure = 64; 401472b9844SJames Wright _pOrder = 3; 402472b9844SJames Wright break; 403472b9844SJames Wright case CGNS_ENUMV(HEXA_125): 404472b9844SJames Wright _numClosure = 125; 405472b9844SJames Wright _pOrder = 4; 406472b9844SJames Wright break; 407472b9844SJames Wright case CGNS_ENUMV(MIXED): 408472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid CGNS ElementType_t: MIXED"); 409472b9844SJames Wright default: 410472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unsupported or Invalid cell type %d", (int)cellType); 411472b9844SJames Wright } 412472b9844SJames Wright if (numClosure) *numClosure = _numClosure; 413472b9844SJames Wright if (pOrder) *pOrder = _pOrder; 414472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 415472b9844SJames Wright } 416472b9844SJames Wright 417472b9844SJames Wright static PetscErrorCode PetscCGNSDataType(PetscDataType pd, CGNS_ENUMT(DataType_t) * cd) 418472b9844SJames Wright { 419472b9844SJames Wright PetscFunctionBegin; 420472b9844SJames Wright switch (pd) { 421472b9844SJames Wright case PETSC_FLOAT: 422472b9844SJames Wright *cd = CGNS_ENUMV(RealSingle); 423472b9844SJames Wright break; 424472b9844SJames Wright case PETSC_DOUBLE: 425472b9844SJames Wright *cd = CGNS_ENUMV(RealDouble); 426472b9844SJames Wright break; 427472b9844SJames Wright case PETSC_COMPLEX: 428472b9844SJames Wright *cd = CGNS_ENUMV(ComplexDouble); 429472b9844SJames Wright break; 430472b9844SJames Wright default: 431472b9844SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Data type %s", PetscDataTypes[pd]); 432472b9844SJames Wright } 433472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 434472b9844SJames Wright } 435472b9844SJames Wright 436472b9844SJames Wright PetscErrorCode DMPlexCreateCGNSFromFile_Internal(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 437472b9844SJames Wright { 438472b9844SJames Wright int cgid = -1; 439472b9844SJames Wright PetscBool use_parallel_viewer = PETSC_FALSE; 440472b9844SJames Wright 441472b9844SJames Wright PetscFunctionBegin; 442472b9844SJames Wright PetscAssertPointer(filename, 2); 4435582b114SJames Wright PetscCall(PetscViewerCGNSRegisterLogEvents_Internal()); 444472b9844SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 445472b9844SJames Wright 446472b9844SJames Wright if (use_parallel_viewer) { 447472b9844SJames Wright PetscCallCGNS(cgp_mpi_comm(comm)); 4484c155f28SJames Wright PetscCallCGNSOpen(cgp_open(filename, CG_MODE_READ, &cgid), 0, 0); 449472b9844SJames Wright PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cgp_open(\"%s\",...) did not return a valid file ID", filename); 450472b9844SJames Wright PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 4514c155f28SJames Wright PetscCallCGNSClose(cgp_close(cgid), 0, 0); 452472b9844SJames Wright } else { 4534c155f28SJames Wright PetscCallCGNSOpen(cg_open(filename, CG_MODE_READ, &cgid), 0, 0); 454472b9844SJames Wright PetscCheck(cgid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 455472b9844SJames Wright PetscCall(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 4564c155f28SJames Wright PetscCallCGNSClose(cg_close(cgid), 0, 0); 457472b9844SJames Wright } 458472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 459472b9844SJames Wright } 460472b9844SJames Wright 461472b9844SJames Wright PetscErrorCode DMPlexCreateCGNS_Internal_Serial(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 462472b9844SJames Wright { 463472b9844SJames Wright PetscMPIInt num_proc, rank; 464472b9844SJames Wright DM cdm; 465472b9844SJames Wright DMLabel label; 466472b9844SJames Wright PetscSection coordSection; 467472b9844SJames Wright Vec coordinates; 468472b9844SJames Wright PetscScalar *coords; 469472b9844SJames Wright PetscInt *cellStart, *vertStart, v; 470472b9844SJames Wright PetscInt labelIdRange[2], labelId; 471472b9844SJames Wright /* Read from file */ 472472b9844SJames Wright char basename[CGIO_MAX_NAME_LENGTH + 1]; 473472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 474472b9844SJames Wright int dim = 0, physDim = 0, coordDim = 0, numVertices = 0, numCells = 0; 475472b9844SJames Wright int nzones = 0; 476472b9844SJames Wright const int B = 1; // Only support single base 477472b9844SJames Wright 478472b9844SJames Wright PetscFunctionBegin; 479472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 480472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 481472b9844SJames Wright PetscCall(DMCreate(comm, dm)); 482472b9844SJames Wright PetscCall(DMSetType(*dm, DMPLEX)); 483472b9844SJames Wright 484472b9844SJames Wright /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 485472b9844SJames Wright if (rank == 0) { 486472b9844SJames Wright int nbases, z; 487472b9844SJames Wright 4884c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 489472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 4904c155f28SJames Wright PetscCallCGNSRead(cg_base_read(cgid, B, basename, &dim, &physDim), *dm, 0); 4914c155f28SJames Wright PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), *dm, 0); 492472b9844SJames Wright PetscCall(PetscCalloc2(nzones + 1, &cellStart, nzones + 1, &vertStart)); 493472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 494472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 495472b9844SJames Wright 4964c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 497472b9844SJames Wright numVertices += sizes[0]; 498472b9844SJames Wright numCells += sizes[1]; 499472b9844SJames Wright cellStart[z] += sizes[1] + cellStart[z - 1]; 500472b9844SJames Wright vertStart[z] += sizes[0] + vertStart[z - 1]; 501472b9844SJames Wright } 502472b9844SJames Wright for (z = 1; z <= nzones; ++z) vertStart[z] += numCells; 503472b9844SJames Wright coordDim = dim; 504472b9844SJames Wright } 505472b9844SJames Wright PetscCallMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH + 1, MPI_CHAR, 0, comm)); 506472b9844SJames Wright PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 507472b9844SJames Wright PetscCallMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm)); 508472b9844SJames Wright PetscCallMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm)); 509472b9844SJames Wright 510472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 511472b9844SJames Wright PetscCall(DMSetDimension(*dm, dim)); 512472b9844SJames Wright PetscCall(DMCreateLabel(*dm, "celltype")); 513472b9844SJames Wright PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 514472b9844SJames Wright 515472b9844SJames Wright /* Read zone information */ 516472b9844SJames Wright if (rank == 0) { 517472b9844SJames Wright int z, c, c_loc; 518472b9844SJames Wright 519472b9844SJames Wright /* Read the cell set connectivity table and build mesh topology 520472b9844SJames Wright CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 521472b9844SJames Wright /* First set sizes */ 522472b9844SJames Wright for (z = 1, c = 0; z <= nzones; ++z) { 523472b9844SJames Wright CGNS_ENUMT(ZoneType_t) zonetype; 524472b9844SJames Wright int nsections; 525472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 526472b9844SJames Wright cgsize_t start, end; 527472b9844SJames Wright int nbndry, parentFlag; 528472b9844SJames Wright PetscInt numCorners, pOrder; 529472b9844SJames Wright DMPolytopeType ctype; 530472b9844SJames Wright const int S = 1; // Only support single section 531472b9844SJames Wright 5324c155f28SJames Wright PetscCallCGNSRead(cg_zone_type(cgid, B, z, &zonetype), *dm, 0); 533472b9844SJames Wright PetscCheck(zonetype != CGNS_ENUMV(Structured), PETSC_COMM_SELF, PETSC_ERR_LIB, "Can only handle Unstructured zones for CGNS"); 5344c155f28SJames Wright PetscCallCGNSRead(cg_nsections(cgid, B, z, &nsections), *dm, 0); 535472b9844SJames Wright PetscCheck(nsections <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single section, not %d", nsections); 5364c155f28SJames Wright PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 537472b9844SJames Wright if (cellType == CGNS_ENUMV(MIXED)) { 538472b9844SJames Wright cgsize_t elementDataSize, *elements; 539472b9844SJames Wright PetscInt off; 540472b9844SJames Wright 5414c155f28SJames Wright PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 542472b9844SJames Wright PetscCall(PetscMalloc1(elementDataSize, &elements)); 5434c155f28SJames Wright PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 544472b9844SJames Wright for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 545ca328833SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[off], &ctype, &numCorners, NULL)); 546ca328833SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[off], NULL, &pOrder)); 547472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 548472b9844SJames Wright PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 549472b9844SJames Wright PetscCall(DMPlexSetCellType(*dm, c, ctype)); 550472b9844SJames Wright off += numCorners + 1; 551472b9844SJames Wright } 552472b9844SJames Wright PetscCall(PetscFree(elements)); 553472b9844SJames Wright } else { 554472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &ctype, &numCorners, NULL)); 555472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 556472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 557472b9844SJames Wright for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 558472b9844SJames Wright PetscCall(DMPlexSetConeSize(*dm, c, numCorners)); 559472b9844SJames Wright PetscCall(DMPlexSetCellType(*dm, c, ctype)); 560472b9844SJames Wright } 561472b9844SJames Wright } 562472b9844SJames Wright } 563472b9844SJames Wright for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 564472b9844SJames Wright } 565472b9844SJames Wright 566472b9844SJames Wright PetscCall(DMSetUp(*dm)); 567472b9844SJames Wright 568472b9844SJames Wright PetscCall(DMCreateLabel(*dm, "zone")); 569472b9844SJames Wright if (rank == 0) { 570472b9844SJames Wright int z, c, c_loc, v_loc; 571472b9844SJames Wright 572472b9844SJames Wright PetscCall(DMGetLabel(*dm, "zone", &label)); 573472b9844SJames Wright for (z = 1, c = 0; z <= nzones; ++z) { 574472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 575472b9844SJames Wright cgsize_t elementDataSize, *elements, start, end; 576472b9844SJames Wright int nbndry, parentFlag; 577472b9844SJames Wright PetscInt *cone, numc, numCorners, maxCorners = 27, pOrder; 578472b9844SJames Wright const int S = 1; // Only support single section 579472b9844SJames Wright 5804c155f28SJames Wright PetscCallCGNSRead(cg_section_read(cgid, B, z, S, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 581472b9844SJames Wright numc = end - start; 5824c155f28SJames Wright PetscCallCGNSRead(cg_ElementDataSize(cgid, B, z, S, &elementDataSize), *dm, 0); 583472b9844SJames Wright PetscCall(PetscMalloc2(elementDataSize, &elements, maxCorners, &cone)); 5844c155f28SJames Wright PetscCallCGNSReadData(cg_poly_elements_read(cgid, B, z, S, elements, NULL, NULL), *dm, 0); 585472b9844SJames Wright if (cellType == CGNS_ENUMV(MIXED)) { 586472b9844SJames Wright /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 587472b9844SJames Wright for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 588ca328833SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &numCorners, NULL)); 589ca328833SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo((CGNS_ENUMT(ElementType_t))elements[v], NULL, &pOrder)); 590472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 591472b9844SJames Wright ++v; 592472b9844SJames Wright for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 593472b9844SJames Wright PetscCall(DMPlexReorderCell(*dm, c, cone)); 594472b9844SJames Wright PetscCall(DMPlexSetCone(*dm, c, cone)); 595472b9844SJames Wright PetscCall(DMLabelSetValue(label, c, z)); 596472b9844SJames Wright } 597472b9844SJames Wright } else { 598472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, NULL, &numCorners, NULL)); 599472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, NULL, &pOrder)); 600472b9844SJames Wright PetscCheck(pOrder == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Serial CGNS reader only supports first order elements, not %" PetscInt_FMT " order", pOrder); 601472b9844SJames Wright /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 602472b9844SJames Wright for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 603472b9844SJames Wright for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) cone[v_loc] = elements[v] + numCells - 1; 604472b9844SJames Wright PetscCall(DMPlexReorderCell(*dm, c, cone)); 605472b9844SJames Wright PetscCall(DMPlexSetCone(*dm, c, cone)); 606472b9844SJames Wright PetscCall(DMLabelSetValue(label, c, z)); 607472b9844SJames Wright } 608472b9844SJames Wright } 609472b9844SJames Wright PetscCall(PetscFree2(elements, cone)); 610472b9844SJames Wright } 611472b9844SJames Wright } 612472b9844SJames Wright 613472b9844SJames Wright PetscCall(DMPlexSymmetrize(*dm)); 614472b9844SJames Wright PetscCall(DMPlexStratify(*dm)); 61501432a30SJames Wright if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 616472b9844SJames Wright 617472b9844SJames Wright /* Read coordinates */ 618472b9844SJames Wright PetscCall(DMSetCoordinateDim(*dm, coordDim)); 619472b9844SJames Wright PetscCall(DMGetCoordinateDM(*dm, &cdm)); 620472b9844SJames Wright PetscCall(DMGetLocalSection(cdm, &coordSection)); 621472b9844SJames Wright PetscCall(PetscSectionSetNumFields(coordSection, 1)); 622472b9844SJames Wright PetscCall(PetscSectionSetFieldComponents(coordSection, 0, coordDim)); 623472b9844SJames Wright PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 624472b9844SJames Wright for (v = numCells; v < numCells + numVertices; ++v) { 625472b9844SJames Wright PetscCall(PetscSectionSetDof(coordSection, v, dim)); 626472b9844SJames Wright PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, coordDim)); 627472b9844SJames Wright } 628472b9844SJames Wright PetscCall(PetscSectionSetUp(coordSection)); 629472b9844SJames Wright 630472b9844SJames Wright PetscCall(DMCreateLocalVector(cdm, &coordinates)); 631472b9844SJames Wright PetscCall(VecGetArray(coordinates, &coords)); 632472b9844SJames Wright if (rank == 0) { 633472b9844SJames Wright PetscInt off = 0; 634472b9844SJames Wright float *x[3]; 635472b9844SJames Wright int z, d; 636472b9844SJames Wright 637472b9844SJames Wright PetscCall(PetscMalloc3(numVertices, &x[0], numVertices, &x[1], numVertices, &x[2])); 638472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 639472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 640472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 641472b9844SJames Wright cgsize_t range_min[3] = {1, 1, 1}; 642472b9844SJames Wright cgsize_t range_max[3] = {1, 1, 1}; 643472b9844SJames Wright int ngrids, ncoords; 644472b9844SJames Wright 6454c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), *dm, 0); 646472b9844SJames Wright range_max[0] = sizes[0]; 6474c155f28SJames Wright PetscCallCGNSRead(cg_ngrids(cgid, B, z, &ngrids), *dm, 0); 648472b9844SJames Wright PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 6494c155f28SJames Wright PetscCallCGNSRead(cg_ncoords(cgid, B, z, &ncoords), *dm, 0); 650472b9844SJames Wright PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 651472b9844SJames Wright for (d = 0; d < coordDim; ++d) { 6524c155f28SJames Wright PetscCallCGNSRead(cg_coord_info(cgid, B, z, 1 + d, &datatype, buffer), *dm, 0); 6534c155f28SJames Wright PetscCallCGNSReadData(cg_coord_read(cgid, B, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]), *dm, 0); 654472b9844SJames Wright } 655472b9844SJames Wright if (coordDim >= 1) { 656472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 0] = x[0][v]; 657472b9844SJames Wright } 658472b9844SJames Wright if (coordDim >= 2) { 659472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 1] = x[1][v]; 660472b9844SJames Wright } 661472b9844SJames Wright if (coordDim >= 3) { 662472b9844SJames Wright for (v = 0; v < sizes[0]; ++v) coords[(v + off) * coordDim + 2] = x[2][v]; 663472b9844SJames Wright } 664472b9844SJames Wright off += sizes[0]; 665472b9844SJames Wright } 666472b9844SJames Wright PetscCall(PetscFree3(x[0], x[1], x[2])); 667472b9844SJames Wright } 668472b9844SJames Wright PetscCall(VecRestoreArray(coordinates, &coords)); 669472b9844SJames Wright 670472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 671472b9844SJames Wright PetscCall(VecSetBlockSize(coordinates, coordDim)); 672472b9844SJames Wright PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 673472b9844SJames Wright PetscCall(VecDestroy(&coordinates)); 674472b9844SJames Wright 675472b9844SJames Wright /* Read boundary conditions */ 676472b9844SJames Wright PetscCall(DMGetNumLabels(*dm, &labelIdRange[0])); 677472b9844SJames Wright if (rank == 0) { 678472b9844SJames Wright CGNS_ENUMT(BCType_t) bctype; 679472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 680472b9844SJames Wright CGNS_ENUMT(PointSetType_t) pointtype; 681472b9844SJames Wright cgsize_t *points; 682472b9844SJames Wright PetscReal *normals; 683472b9844SJames Wright int normal[3]; 684472b9844SJames Wright char *bcname = buffer; 685472b9844SJames Wright cgsize_t npoints, nnormals; 686472b9844SJames Wright int z, nbc, bc, c, ndatasets; 687472b9844SJames Wright 688472b9844SJames Wright for (z = 1; z <= nzones; ++z) { 6894c155f28SJames Wright PetscCallCGNSRead(cg_nbocos(cgid, B, z, &nbc), *dm, 0); 690472b9844SJames Wright for (bc = 1; bc <= nbc; ++bc) { 6914c155f28SJames Wright PetscCallCGNSRead(cg_boco_info(cgid, B, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets), *dm, 0); 692472b9844SJames Wright PetscCall(DMCreateLabel(*dm, bcname)); 693472b9844SJames Wright PetscCall(DMGetLabel(*dm, bcname, &label)); 694472b9844SJames Wright PetscCall(PetscMalloc2(npoints, &points, nnormals, &normals)); 6954c155f28SJames Wright PetscCallCGNSReadData(cg_boco_read(cgid, B, z, bc, points, (void *)normals), *dm, 0); 696472b9844SJames Wright if (pointtype == CGNS_ENUMV(ElementRange)) { 697472b9844SJames Wright // Range of cells: assuming half-open interval 698472b9844SJames Wright for (c = points[0]; c < points[1]; ++c) PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 699472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(ElementList)) { 700472b9844SJames Wright // List of cells 701472b9844SJames Wright for (c = 0; c < npoints; ++c) PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 702472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(PointRange)) { 703472b9844SJames Wright CGNS_ENUMT(GridLocation_t) gridloc; 704472b9844SJames Wright 705472b9844SJames Wright // List of points: 706472b9844SJames Wright PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 7074c155f28SJames Wright PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 708472b9844SJames Wright // Range of points: assuming half-open interval 709472b9844SJames Wright for (c = points[0]; c < points[1]; ++c) { 710472b9844SJames Wright if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, c - vertStart[z - 1], 1)); 711472b9844SJames Wright else PetscCall(DMLabelSetValue(label, c - cellStart[z - 1], 1)); 712472b9844SJames Wright } 713472b9844SJames Wright } else if (pointtype == CGNS_ENUMV(PointList)) { 714472b9844SJames Wright CGNS_ENUMT(GridLocation_t) gridloc; 715472b9844SJames Wright 716472b9844SJames Wright // List of points: 717472b9844SJames Wright PetscCallCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 7184c155f28SJames Wright PetscCallCGNSRead(cg_gridlocation_read(&gridloc), *dm, 0); 719472b9844SJames Wright for (c = 0; c < npoints; ++c) { 720472b9844SJames Wright if (gridloc == CGNS_ENUMV(Vertex)) PetscCall(DMLabelSetValue(label, points[c] - vertStart[z - 1], 1)); 721472b9844SJames Wright else PetscCall(DMLabelSetValue(label, points[c] - cellStart[z - 1], 1)); 722472b9844SJames Wright } 723472b9844SJames Wright } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int)pointtype); 724472b9844SJames Wright PetscCall(PetscFree2(points, normals)); 725472b9844SJames Wright } 726472b9844SJames Wright } 727472b9844SJames Wright PetscCall(PetscFree2(cellStart, vertStart)); 728472b9844SJames Wright } 729472b9844SJames Wright PetscCall(DMGetNumLabels(*dm, &labelIdRange[1])); 730472b9844SJames Wright PetscCallMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm)); 731472b9844SJames Wright 732472b9844SJames Wright /* Create BC labels at all processes */ 733472b9844SJames Wright for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 734472b9844SJames Wright char *labelName = buffer; 735472b9844SJames Wright size_t len = sizeof(buffer); 736472b9844SJames Wright const char *locName; 737472b9844SJames Wright 738472b9844SJames Wright if (rank == 0) { 739472b9844SJames Wright PetscCall(DMGetLabelByNum(*dm, labelId, &label)); 740472b9844SJames Wright PetscCall(PetscObjectGetName((PetscObject)label, &locName)); 741472b9844SJames Wright PetscCall(PetscStrncpy(labelName, locName, len)); 742472b9844SJames Wright } 743472b9844SJames Wright PetscCallMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm)); 744472b9844SJames Wright PetscCallMPI(DMCreateLabel(*dm, labelName)); 745472b9844SJames Wright } 746472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 747472b9844SJames Wright } 748472b9844SJames Wright 749472b9844SJames Wright PetscErrorCode DMPlexCreateCGNS_Internal_Parallel(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 750472b9844SJames Wright { 751472b9844SJames Wright PetscMPIInt num_proc, rank; 752472b9844SJames Wright /* Read from file */ 753472b9844SJames Wright char basename[CGIO_MAX_NAME_LENGTH + 1]; 754472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 755472b9844SJames Wright int dim = 0, physDim = 0, coordDim = 0; 756472b9844SJames Wright PetscInt NVertices = 0, NCells = 0; 757472b9844SJames Wright int nzones = 0, nbases; 758b30057acSJames Wright int zone = 1; // Only supports single zone files 759b30057acSJames Wright int base = 1; // Only supports single base 760472b9844SJames Wright 761472b9844SJames Wright PetscFunctionBegin; 762472b9844SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 763472b9844SJames Wright PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 764472b9844SJames Wright PetscCall(DMCreate(comm, dm)); 765472b9844SJames Wright PetscCall(DMSetType(*dm, DMPLEX)); 766472b9844SJames Wright 7674c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), *dm, 0); 768472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 769f605775fSPierre Jolivet // From the CGNS web page cell_dim phys_dim (embedding space in PETSc) CGNS defines as length of spatial vectors/components) 770b30057acSJames Wright PetscCallCGNSRead(cg_base_read(cgid, base, basename, &dim, &physDim), *dm, 0); 771b30057acSJames Wright PetscCallCGNSRead(cg_nzones(cgid, base, &nzones), *dm, 0); 772472b9844SJames Wright PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "Parallel reader limited to one zone, not %d", nzones); 773472b9844SJames Wright { 774472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 775472b9844SJames Wright 776b30057acSJames Wright PetscCallCGNSRead(cg_zone_read(cgid, base, zone, buffer, sizes), *dm, 0); 777472b9844SJames Wright NVertices = sizes[0]; 778472b9844SJames Wright NCells = sizes[1]; 779472b9844SJames Wright } 780472b9844SJames Wright 781472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm, basename)); 782472b9844SJames Wright PetscCall(DMSetDimension(*dm, dim)); 783472b9844SJames Wright coordDim = dim; 784472b9844SJames Wright 785472b9844SJames 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 786472b9844SJames Wright // call to get this right but continuing for now with single section, single topology, one zone. 787472b9844SJames Wright // establish element ranges for my rank 788472b9844SJames Wright PetscInt mystarte, myende, mystartv, myendv, myownede, myownedv; 789472b9844SJames Wright PetscLayout elem_map, vtx_map; 790472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NCells, 1, &elem_map)); 791472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 792472b9844SJames Wright PetscCall(PetscLayoutGetRange(elem_map, &mystarte, &myende)); 793472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(elem_map, &myownede)); 794472b9844SJames Wright PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 795472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 796472b9844SJames Wright 797472b9844SJames Wright // -- Build Plex in parallel 798f1a88294SStefano Zampini DMPolytopeType dm_cell_type = DM_POLYTOPE_UNKNOWN; 799472b9844SJames Wright PetscInt pOrder = 1, numClosure = -1; 800472b9844SJames Wright cgsize_t *elements; 801472b9844SJames Wright { 802472b9844SJames Wright int nsections; 803472b9844SJames Wright PetscInt *elementsQ1, numCorners = -1; 804472b9844SJames Wright const int *perm; 805472b9844SJames Wright cgsize_t start, end; // Throwaway 806472b9844SJames Wright 807b30057acSJames Wright cg_nsections(cgid, base, zone, &nsections); 808472b9844SJames Wright // Read element connectivity 809472b9844SJames Wright for (int index_sect = 1; index_sect <= nsections; index_sect++) { 810472b9844SJames Wright int nbndry, parentFlag; 811472b9844SJames Wright PetscInt cell_dim; 812472b9844SJames Wright CGNS_ENUMT(ElementType_t) cellType; 813472b9844SJames Wright 814b30057acSJames Wright PetscCallCGNSRead(cg_section_read(cgid, base, zone, index_sect, buffer, &cellType, &start, &end, &nbndry, &parentFlag), *dm, 0); 815472b9844SJames Wright 816472b9844SJames Wright PetscCall(CGNSElementTypeGetTopologyInfo(cellType, &dm_cell_type, &numCorners, &cell_dim)); 817472b9844SJames Wright // Skip over element that are not max dimension (ie. boundary elements) 818472b9844SJames Wright if (cell_dim != dim) continue; 819472b9844SJames Wright PetscCall(CGNSElementTypeGetDiscretizationInfo(cellType, &numClosure, &pOrder)); 820472b9844SJames Wright PetscCall(PetscMalloc1(myownede * numClosure, &elements)); 821b30057acSJames Wright PetscCallCGNSReadData(cgp_elements_read_data(cgid, base, zone, index_sect, mystarte + 1, myende, elements), *dm, 0); 822472b9844SJames Wright for (PetscInt v = 0; v < myownede * numClosure; ++v) elements[v] -= 1; // 0 based 823472b9844SJames Wright break; 824472b9844SJames Wright } 825472b9844SJames Wright 826472b9844SJames Wright // Create corners-only connectivity 827472b9844SJames Wright PetscCall(PetscMalloc1(myownede * numCorners, &elementsQ1)); 828472b9844SJames Wright PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numCorners, NULL, &perm)); 829472b9844SJames Wright for (PetscInt e = 0; e < myownede; ++e) { 830472b9844SJames Wright for (PetscInt v = 0; v < numCorners; ++v) elementsQ1[e * numCorners + perm[v]] = elements[e * numClosure + v]; 831472b9844SJames Wright } 832472b9844SJames Wright 833472b9844SJames Wright // Build cell-vertex Plex 834472b9844SJames Wright PetscCall(DMPlexBuildFromCellListParallel(*dm, myownede, myownedv, NVertices, numCorners, elementsQ1, NULL, NULL)); 835472b9844SJames Wright PetscCall(DMViewFromOptions(*dm, NULL, "-corner_dm_view")); 836472b9844SJames Wright PetscCall(PetscFree(elementsQ1)); 837472b9844SJames Wright } 838472b9844SJames Wright 83901432a30SJames Wright if (interpolate) PetscCall(DMPlexInterpolateInPlace_Internal(*dm)); 840472b9844SJames Wright 841472b9844SJames Wright // -- Create SF for naive nodal-data read to elements 842472b9844SJames Wright PetscSF plex_to_cgns_sf; 843472b9844SJames Wright { 844472b9844SJames Wright PetscInt nleaves, num_comp; 845472b9844SJames Wright PetscInt *leaf, num_leaves = 0; 846472b9844SJames Wright PetscInt cStart, cEnd; 847472b9844SJames Wright const int *perm; 848472b9844SJames Wright PetscSF cgns_to_local_sf; 849472b9844SJames Wright PetscSection local_section; 850472b9844SJames Wright PetscFE fe; 851472b9844SJames Wright 852472b9844SJames Wright // sfNatural requires PetscSection to handle DMDistribute, so we use PetscFE to define the section 853472b9844SJames Wright // Use number of components = 1 to work with just the nodes themselves 854472b9844SJames Wright PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, 1, dm_cell_type, pOrder, PETSC_DETERMINE, &fe)); 855472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)fe, "FE for sfNatural")); 856472b9844SJames Wright PetscCall(DMAddField(*dm, NULL, (PetscObject)fe)); 857472b9844SJames Wright PetscCall(DMCreateDS(*dm)); 858472b9844SJames Wright PetscCall(PetscFEDestroy(&fe)); 859472b9844SJames Wright 860472b9844SJames Wright PetscCall(DMGetLocalSection(*dm, &local_section)); 861472b9844SJames Wright PetscCall(PetscSectionViewFromOptions(local_section, NULL, "-fe_natural_section_view")); 862472b9844SJames Wright PetscCall(PetscSectionGetFieldComponents(local_section, 0, &num_comp)); 863472b9844SJames Wright PetscCall(PetscSectionGetStorageSize(local_section, &nleaves)); 864472b9844SJames Wright nleaves /= num_comp; 865472b9844SJames Wright PetscCall(PetscMalloc1(nleaves, &leaf)); 866472b9844SJames Wright for (PetscInt i = 0; i < nleaves; i++) leaf[i] = -1; 867472b9844SJames Wright 868472b9844SJames Wright // Get the permutation from CGNS closure numbering to PLEX closure numbering 869472b9844SJames Wright PetscCall(DMPlexCGNSGetPermutation_Internal(dm_cell_type, numClosure, NULL, &perm)); 870472b9844SJames Wright PetscCall(DMPlexGetHeightStratum(*dm, 0, &cStart, &cEnd)); 871472b9844SJames Wright for (PetscInt cell = cStart; cell < cEnd; ++cell) { 872472b9844SJames Wright PetscInt num_closure_dof, *closure_idx = NULL; 873472b9844SJames Wright 874472b9844SJames Wright PetscCall(DMPlexGetClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 875472b9844SJames Wright PetscAssert(numClosure * num_comp == num_closure_dof, comm, PETSC_ERR_PLIB, "Closure dof size does not match polytope"); 876472b9844SJames Wright for (PetscInt i = 0; i < numClosure; i++) { 877472b9844SJames Wright PetscInt li = closure_idx[perm[i] * num_comp] / num_comp; 878472b9844SJames Wright if (li < 0) continue; 879472b9844SJames Wright 880472b9844SJames Wright PetscInt cgns_idx = elements[cell * numClosure + i]; 881472b9844SJames Wright if (leaf[li] == -1) { 882472b9844SJames Wright leaf[li] = cgns_idx; 883472b9844SJames Wright num_leaves++; 884472b9844SJames Wright } else PetscAssert(leaf[li] == cgns_idx, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf does not match previously set"); 885472b9844SJames Wright } 886472b9844SJames Wright PetscCall(DMPlexRestoreClosureIndices(*dm, local_section, local_section, cell, PETSC_FALSE, &num_closure_dof, &closure_idx, NULL, NULL)); 887472b9844SJames Wright } 888472b9844SJames Wright PetscAssert(num_leaves == nleaves, PETSC_COMM_SELF, PETSC_ERR_PLIB, "leaf count in closure does not match nleaves"); 889472b9844SJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)*dm), &cgns_to_local_sf)); 890472b9844SJames Wright PetscCall(PetscSFSetGraphLayout(cgns_to_local_sf, vtx_map, nleaves, NULL, PETSC_USE_POINTER, leaf)); 891472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)cgns_to_local_sf, "CGNS to Plex SF")); 892472b9844SJames Wright PetscCall(PetscSFViewFromOptions(cgns_to_local_sf, NULL, "-CGNStoPlex_sf_view")); 893472b9844SJames Wright PetscCall(PetscFree(leaf)); 894472b9844SJames Wright PetscCall(PetscFree(elements)); 895472b9844SJames Wright 896472b9844SJames Wright { // Convert cgns_to_local to global_to_cgns 897472b9844SJames Wright PetscSF sectionsf, cgns_to_global_sf; 898472b9844SJames Wright 899472b9844SJames Wright PetscCall(DMGetSectionSF(*dm, §ionsf)); 900472b9844SJames Wright PetscCall(PetscSFComposeInverse(cgns_to_local_sf, sectionsf, &cgns_to_global_sf)); 901472b9844SJames Wright PetscCall(PetscSFDestroy(&cgns_to_local_sf)); 902472b9844SJames Wright PetscCall(PetscSFCreateInverseSF(cgns_to_global_sf, &plex_to_cgns_sf)); 903472b9844SJames Wright PetscCall(PetscObjectSetName((PetscObject)plex_to_cgns_sf, "Global Plex to CGNS")); 904472b9844SJames Wright PetscCall(PetscSFDestroy(&cgns_to_global_sf)); 905472b9844SJames Wright } 906472b9844SJames Wright } 907472b9844SJames Wright 908472b9844SJames Wright { // -- Set coordinates for DM 909472b9844SJames Wright PetscScalar *coords; 910472b9844SJames Wright float *x[3]; 911472b9844SJames Wright double *xd[3]; 912472b9844SJames Wright PetscBool read_with_double; 913472b9844SJames Wright PetscFE cfe; 914472b9844SJames Wright 915472b9844SJames Wright // Setup coordinate space first. Use pOrder here for isoparametric; revisit with CPEX-0045 High Order. 916472b9844SJames Wright PetscCall(PetscFECreateLagrangeByCell(PETSC_COMM_SELF, dim, coordDim, dm_cell_type, pOrder, PETSC_DETERMINE, &cfe)); 9174c712d99Sksagiyam PetscCall(DMSetCoordinateDisc(*dm, cfe, PETSC_FALSE, PETSC_FALSE)); 918472b9844SJames Wright PetscCall(PetscFEDestroy(&cfe)); 919472b9844SJames Wright 920472b9844SJames Wright { // Determine if coords are written in single or double precision 921472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 922472b9844SJames Wright 923b30057acSJames Wright PetscCallCGNSRead(cg_coord_info(cgid, base, zone, 1, &datatype, buffer), *dm, 0); 9245582b114SJames Wright read_with_double = datatype == CGNS_ENUMV(RealDouble) ? PETSC_TRUE : PETSC_FALSE; 925472b9844SJames Wright } 926472b9844SJames Wright 927472b9844SJames Wright // Read coords from file and set into component-major ordering 928472b9844SJames Wright if (read_with_double) PetscCall(PetscMalloc3(myownedv, &xd[0], myownedv, &xd[1], myownedv, &xd[2])); 929472b9844SJames Wright else PetscCall(PetscMalloc3(myownedv, &x[0], myownedv, &x[1], myownedv, &x[2])); 930472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * coordDim, &coords)); 931472b9844SJames Wright { 932472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 933472b9844SJames Wright cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 934472b9844SJames Wright cgsize_t range_max[3] = {myendv, 1, 1}; 935472b9844SJames Wright int ngrids, ncoords; 936472b9844SJames Wright 937b30057acSJames Wright PetscCallCGNSRead(cg_zone_read(cgid, base, zone, buffer, sizes), *dm, 0); 938b30057acSJames Wright PetscCallCGNSRead(cg_ngrids(cgid, base, zone, &ngrids), *dm, 0); 939472b9844SJames Wright PetscCheck(ngrids <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single grid, not %d", ngrids); 940b30057acSJames Wright PetscCallCGNSRead(cg_ncoords(cgid, base, zone, &ncoords), *dm, 0); 941472b9844SJames Wright PetscCheck(ncoords == coordDim, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a coordinate array for each dimension, not %d", ncoords); 942472b9844SJames Wright if (read_with_double) { 943b30057acSJames Wright for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, base, zone, (d + 1), range_min, range_max, xd[d]), *dm, 0); 944472b9844SJames Wright if (coordDim >= 1) { 945472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = xd[0][v]; 946472b9844SJames Wright } 947472b9844SJames Wright if (coordDim >= 2) { 948472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = xd[1][v]; 949472b9844SJames Wright } 950472b9844SJames Wright if (coordDim >= 3) { 951472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = xd[2][v]; 952472b9844SJames Wright } 953472b9844SJames Wright } else { 954b30057acSJames Wright for (int d = 0; d < coordDim; ++d) PetscCallCGNSReadData(cgp_coord_read_data(cgid, 1, zone, (d + 1), range_min, range_max, x[d]), *dm, 0); 955472b9844SJames Wright if (coordDim >= 1) { 956472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 0] = x[0][v]; 957472b9844SJames Wright } 958472b9844SJames Wright if (coordDim >= 2) { 959472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 1] = x[1][v]; 960472b9844SJames Wright } 961472b9844SJames Wright if (coordDim >= 3) { 962472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) coords[v * coordDim + 2] = x[2][v]; 963472b9844SJames Wright } 964472b9844SJames Wright } 965472b9844SJames Wright } 966472b9844SJames Wright if (read_with_double) PetscCall(PetscFree3(xd[0], xd[1], xd[2])); 967472b9844SJames Wright else PetscCall(PetscFree3(x[0], x[1], x[2])); 968472b9844SJames Wright 969472b9844SJames Wright { // Reduce CGNS-ordered coordinate nodes to Plex ordering, and set DM's coordinates 970472b9844SJames Wright Vec coord_global; 971472b9844SJames Wright MPI_Datatype unit; 972472b9844SJames Wright PetscScalar *coord_global_array; 973472b9844SJames Wright DM cdm; 974472b9844SJames Wright 975472b9844SJames Wright PetscCall(DMGetCoordinateDM(*dm, &cdm)); 976472b9844SJames Wright PetscCall(DMCreateGlobalVector(cdm, &coord_global)); 977472b9844SJames Wright PetscCall(VecGetArrayWrite(coord_global, &coord_global_array)); 978472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(coordDim, MPIU_SCALAR, &unit)); 979472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 980472b9844SJames Wright PetscCall(PetscSFReduceBegin(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 981472b9844SJames Wright PetscCall(PetscSFReduceEnd(plex_to_cgns_sf, unit, coords, coord_global_array, MPI_REPLACE)); 982472b9844SJames Wright PetscCall(VecRestoreArrayWrite(coord_global, &coord_global_array)); 983472b9844SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 984472b9844SJames Wright PetscCall(DMSetCoordinates(*dm, coord_global)); 985472b9844SJames Wright PetscCall(VecDestroy(&coord_global)); 986472b9844SJames Wright } 987472b9844SJames Wright PetscCall(PetscFree(coords)); 988472b9844SJames Wright } 989472b9844SJames Wright 990472b9844SJames Wright // -- Set sfNatural for solution vectors in CGNS file 991472b9844SJames Wright // NOTE: We set sfNatural to be the map between the original CGNS ordering of nodes and the Plex ordering of nodes. 992472b9844SJames Wright PetscCall(PetscSFViewFromOptions(plex_to_cgns_sf, NULL, "-sfNatural_init_view")); 993472b9844SJames Wright PetscCall(DMSetNaturalSF(*dm, plex_to_cgns_sf)); 994472b9844SJames Wright PetscCall(DMSetUseNatural(*dm, PETSC_TRUE)); 995472b9844SJames Wright PetscCall(PetscSFDestroy(&plex_to_cgns_sf)); 996472b9844SJames Wright 997472b9844SJames Wright PetscCall(PetscLayoutDestroy(&elem_map)); 998472b9844SJames Wright PetscCall(PetscLayoutDestroy(&vtx_map)); 9993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10005f34f2dcSJed Brown } 10015f34f2dcSJed Brown 10025f34f2dcSJed Brown // node_l2g must be freed 1003d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateNodeNumbering(DM dm, PetscInt *num_local_nodes, PetscInt *num_global_nodes, PetscInt *nStart, PetscInt *nEnd, const PetscInt **node_l2g) 1004d71ae5a4SJacob Faibussowitsch { 10055f34f2dcSJed Brown PetscSection local_section; 10065f34f2dcSJed Brown PetscSF point_sf; 10075f34f2dcSJed Brown PetscInt pStart, pEnd, spStart, spEnd, *points, nleaves, ncomp, *nodes; 10085f34f2dcSJed Brown PetscMPIInt comm_size; 10095f34f2dcSJed Brown const PetscInt *ilocal, field = 0; 10105f34f2dcSJed Brown 10115f34f2dcSJed Brown PetscFunctionBegin; 10125f34f2dcSJed Brown *num_local_nodes = -1; 10135f34f2dcSJed Brown *num_global_nodes = -1; 10145f34f2dcSJed Brown *nStart = -1; 10155f34f2dcSJed Brown *nEnd = -1; 10165f34f2dcSJed Brown *node_l2g = NULL; 10175f34f2dcSJed Brown 10185f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &comm_size)); 10195f34f2dcSJed Brown PetscCall(DMGetLocalSection(dm, &local_section)); 10205f34f2dcSJed Brown PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 10215f34f2dcSJed Brown PetscCall(PetscSectionGetChart(local_section, &spStart, &spEnd)); 10225f34f2dcSJed Brown PetscCall(DMGetPointSF(dm, &point_sf)); 10235f34f2dcSJed Brown if (comm_size == 1) nleaves = 0; 10245f34f2dcSJed Brown else PetscCall(PetscSFGetGraph(point_sf, NULL, &nleaves, &ilocal, NULL)); 10255f34f2dcSJed Brown PetscCall(PetscSectionGetFieldComponents(local_section, field, &ncomp)); 10265f34f2dcSJed Brown 10275f34f2dcSJed Brown PetscInt local_node = 0, owned_node = 0, owned_start = 0; 10285f34f2dcSJed Brown for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 10295f34f2dcSJed Brown PetscInt dof; 10305f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10315f34f2dcSJed 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); 10325f34f2dcSJed Brown local_node += dof / ncomp; 10335f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 10345f34f2dcSJed Brown leaf++; 10355f34f2dcSJed Brown } else { 10365f34f2dcSJed Brown owned_node += dof / ncomp; 10375f34f2dcSJed Brown } 10385f34f2dcSJed Brown } 10395f34f2dcSJed Brown PetscCallMPI(MPI_Exscan(&owned_node, &owned_start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 10405f34f2dcSJed Brown PetscCall(PetscMalloc1(pEnd - pStart, &points)); 10415f34f2dcSJed Brown owned_node = 0; 10425f34f2dcSJed Brown for (PetscInt p = spStart, leaf = 0; p < spEnd; p++) { 10435f34f2dcSJed Brown if (leaf < nleaves && p == ilocal[leaf]) { // skip points owned by a different process 10445f34f2dcSJed Brown points[p - pStart] = -1; 10455f34f2dcSJed Brown leaf++; 10465f34f2dcSJed Brown continue; 10475f34f2dcSJed Brown } 10485f34f2dcSJed Brown PetscInt dof, offset; 10495f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10505f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 10515f34f2dcSJed 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); 10525f34f2dcSJed Brown points[p - pStart] = owned_start + owned_node; 10535f34f2dcSJed Brown owned_node += dof / ncomp; 10545f34f2dcSJed Brown } 10555f34f2dcSJed Brown if (comm_size > 1) { 10565f34f2dcSJed Brown PetscCall(PetscSFBcastBegin(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 10575f34f2dcSJed Brown PetscCall(PetscSFBcastEnd(point_sf, MPIU_INT, points, points, MPI_REPLACE)); 10585f34f2dcSJed Brown } 10595f34f2dcSJed Brown 10605f34f2dcSJed Brown // Set up global indices for each local node 10615f34f2dcSJed Brown PetscCall(PetscMalloc1(local_node, &nodes)); 10625f34f2dcSJed Brown for (PetscInt p = spStart; p < spEnd; p++) { 10635f34f2dcSJed Brown PetscInt dof, offset; 10645f34f2dcSJed Brown PetscCall(PetscSectionGetFieldDof(local_section, p, field, &dof)); 10655f34f2dcSJed Brown PetscCall(PetscSectionGetFieldOffset(local_section, p, field, &offset)); 10665f34f2dcSJed Brown for (PetscInt n = 0; n < dof / ncomp; n++) nodes[offset / ncomp + n] = points[p - pStart] + n; 10675f34f2dcSJed Brown } 10685f34f2dcSJed Brown PetscCall(PetscFree(points)); 10695f34f2dcSJed Brown *num_local_nodes = local_node; 10705f34f2dcSJed Brown *nStart = owned_start; 10715f34f2dcSJed Brown *nEnd = owned_start + owned_node; 1072462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&owned_node, num_global_nodes, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 10735f34f2dcSJed Brown *node_l2g = nodes; 10743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10755f34f2dcSJed Brown } 10765f34f2dcSJed Brown 1077d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_PlexCGNS(DM dm, PetscViewer viewer) 1078d71ae5a4SJacob Faibussowitsch { 10795f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 108044698e90SJames Wright PetscInt fvGhostStart; 10815f34f2dcSJed Brown PetscInt topo_dim, coord_dim, num_global_elems; 10825f34f2dcSJed Brown PetscInt cStart, cEnd, num_local_nodes, num_global_nodes, nStart, nEnd; 10835f34f2dcSJed Brown const PetscInt *node_l2g; 10845f34f2dcSJed Brown Vec coord; 1085b85bf5edSJed Brown DM colloc_dm, cdm; 10865f34f2dcSJed Brown PetscMPIInt size; 10875f34f2dcSJed Brown const char *dm_name; 10885f34f2dcSJed Brown int base, zone; 10895f34f2dcSJed Brown cgsize_t isize[3]; 10905f34f2dcSJed Brown 10915f34f2dcSJed Brown PetscFunctionBegin; 109225760affSJed Brown if (!cgv->file_num) { 109325760affSJed Brown PetscInt time_step; 109425760affSJed Brown PetscCall(DMGetOutputSequenceNumber(dm, &time_step, NULL)); 109525760affSJed Brown PetscCall(PetscViewerCGNSFileOpen_Internal(viewer, time_step)); 109625760affSJed Brown } 10975f34f2dcSJed Brown PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size)); 10985f34f2dcSJed Brown PetscCall(DMGetDimension(dm, &topo_dim)); 10995f34f2dcSJed Brown PetscCall(DMGetCoordinateDim(dm, &coord_dim)); 11005f34f2dcSJed Brown PetscCall(PetscObjectGetName((PetscObject)dm, &dm_name)); 11014c155f28SJames Wright PetscCallCGNSWrite(cg_base_write(cgv->file_num, dm_name, topo_dim, coord_dim, &base), dm, viewer); 11025f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, NULL)); 11034c155f28SJames Wright PetscCallCGNSWrite(cg_dataclass_write(CGNS_ENUMV(NormalizedByDimensional)), dm, viewer); 11045f34f2dcSJed Brown 1105b85bf5edSJed Brown { 1106b85bf5edSJed Brown PetscFE fe, fe_coord; 11079bb8f83bSJed Brown PetscClassId ds_id; 1108b85bf5edSJed Brown PetscDualSpace dual_space, dual_space_coord; 11092d1fc720SJed Brown PetscInt num_fields, field_order = -1, field_order_coord; 1110b85bf5edSJed Brown PetscBool is_simplex; 11115b8b2c14SJed Brown PetscCall(DMGetNumFields(dm, &num_fields)); 11129bb8f83bSJed Brown if (num_fields > 0) { 11139bb8f83bSJed Brown PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 11149bb8f83bSJed Brown PetscCall(PetscObjectGetClassId((PetscObject)fe, &ds_id)); 11152d1fc720SJed Brown if (ds_id != PETSCFE_CLASSID) { 11162d1fc720SJed Brown fe = NULL; 11172d1fc720SJed Brown if (ds_id == PETSCFV_CLASSID) field_order = -1; // use whatever is present for coords; field will be CellCenter 11182d1fc720SJed Brown else field_order = 1; // assume vertex-based linear elements 11192d1fc720SJed Brown } 11209bb8f83bSJed Brown } else fe = NULL; 1121b85bf5edSJed Brown if (fe) { 1122b85bf5edSJed Brown PetscCall(PetscFEGetDualSpace(fe, &dual_space)); 1123b85bf5edSJed Brown PetscCall(PetscDualSpaceGetOrder(dual_space, &field_order)); 11242d1fc720SJed Brown } 11255f34f2dcSJed Brown PetscCall(DMGetCoordinateDM(dm, &cdm)); 1126b85bf5edSJed Brown PetscCall(DMGetField(cdm, 0, NULL, (PetscObject *)&fe_coord)); 11279bb8f83bSJed Brown { 11289bb8f83bSJed Brown PetscClassId id; 11299bb8f83bSJed Brown PetscCall(PetscObjectGetClassId((PetscObject)fe_coord, &id)); 11309bb8f83bSJed Brown if (id != PETSCFE_CLASSID) fe_coord = NULL; 11319bb8f83bSJed Brown } 1132b85bf5edSJed Brown if (fe_coord) { 1133b85bf5edSJed Brown PetscCall(PetscFEGetDualSpace(fe_coord, &dual_space_coord)); 1134b85bf5edSJed Brown PetscCall(PetscDualSpaceGetOrder(dual_space_coord, &field_order_coord)); 1135b85bf5edSJed Brown } else field_order_coord = 1; 11362d1fc720SJed Brown if (field_order > 0 && field_order != field_order_coord) { 1137b85bf5edSJed Brown PetscInt quadrature_order = field_order; 1138b85bf5edSJed Brown PetscCall(DMClone(dm, &colloc_dm)); 11397d4eb7abSJed Brown { // Inform the new colloc_dm that it is a coordinate DM so isoperiodic affine corrections can be applied 11401fca310dSJames Wright const PetscSF *face_sfs; 11411fca310dSJames Wright PetscInt num_face_sfs; 11421fca310dSJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(dm, &num_face_sfs, &face_sfs)); 11431fca310dSJames Wright PetscCall(DMPlexSetIsoperiodicFaceSF(colloc_dm, num_face_sfs, (PetscSF *)face_sfs)); 11441fca310dSJames Wright if (face_sfs) colloc_dm->periodic.setup = DMPeriodicCoordinateSetUp_Internal; 11457d4eb7abSJed Brown } 1146b85bf5edSJed Brown PetscCall(DMPlexIsSimplex(dm, &is_simplex)); 1147b85bf5edSJed Brown PetscCall(PetscFECreateLagrange(PetscObjectComm((PetscObject)dm), topo_dim, coord_dim, is_simplex, field_order, quadrature_order, &fe)); 11484c712d99Sksagiyam PetscCall(DMSetCoordinateDisc(colloc_dm, fe, PETSC_FALSE, PETSC_TRUE)); 1149b85bf5edSJed Brown PetscCall(PetscFEDestroy(&fe)); 1150b85bf5edSJed Brown } else { 1151b85bf5edSJed Brown PetscCall(PetscObjectReference((PetscObject)dm)); 1152b85bf5edSJed Brown colloc_dm = dm; 1153b85bf5edSJed Brown } 1154b85bf5edSJed Brown } 1155b85bf5edSJed Brown PetscCall(DMGetCoordinateDM(colloc_dm, &cdm)); 11565f34f2dcSJed Brown PetscCall(DMPlexCreateNodeNumbering(cdm, &num_local_nodes, &num_global_nodes, &nStart, &nEnd, &node_l2g)); 1157b85bf5edSJed Brown PetscCall(DMGetCoordinatesLocal(colloc_dm, &coord)); 11585f34f2dcSJed Brown PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 115944698e90SJames Wright PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 116044698e90SJames Wright if (fvGhostStart >= 0) cEnd = fvGhostStart; 11615f34f2dcSJed Brown num_global_elems = cEnd - cStart; 1162462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(MPI_IN_PLACE, &num_global_elems, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)dm))); 11635f34f2dcSJed Brown isize[0] = num_global_nodes; 11645f34f2dcSJed Brown isize[1] = num_global_elems; 11655f34f2dcSJed Brown isize[2] = 0; 11664c155f28SJames Wright PetscCallCGNSWrite(cg_zone_write(cgv->file_num, base, "Zone", isize, CGNS_ENUMV(Unstructured), &zone), dm, viewer); 11675f34f2dcSJed Brown 11685f34f2dcSJed Brown { 11695f34f2dcSJed Brown const PetscScalar *X; 11705f34f2dcSJed Brown PetscScalar *x; 11715f34f2dcSJed Brown int coord_ids[3]; 11725f34f2dcSJed Brown 11735f34f2dcSJed Brown PetscCall(VecGetArrayRead(coord, &X)); 11745f34f2dcSJed Brown for (PetscInt d = 0; d < coord_dim; d++) { 11755f34f2dcSJed Brown const double exponents[] = {0, 1, 0, 0, 0}; 11765f34f2dcSJed Brown char coord_name[64]; 11775f34f2dcSJed Brown PetscCall(PetscSNPrintf(coord_name, sizeof coord_name, "Coordinate%c", 'X' + (int)d)); 11784c155f28SJames Wright PetscCallCGNSWrite(cgp_coord_write(cgv->file_num, base, zone, CGNS_ENUMV(RealDouble), coord_name, &coord_ids[d]), dm, viewer); 11795f34f2dcSJed Brown PetscCallCGNS(cg_goto(cgv->file_num, base, "Zone_t", zone, "GridCoordinates", 0, coord_name, 0, NULL)); 11804c155f28SJames Wright PetscCallCGNSWrite(cg_exponents_write(CGNS_ENUMV(RealDouble), exponents), dm, viewer); 11815f34f2dcSJed Brown } 11825f34f2dcSJed Brown 11835f34f2dcSJed Brown int section; 11845f34f2dcSJed Brown cgsize_t e_owned, e_global, e_start, *conn = NULL; 11855f34f2dcSJed Brown const int *perm; 1186e853fb4cSJames Wright CGNS_ENUMT(ElementType_t) element_type = CGNS_ENUMV(ElementTypeNull); 11875f34f2dcSJed Brown { 11885f34f2dcSJed Brown PetscCall(PetscMalloc1(nEnd - nStart, &x)); 11895f34f2dcSJed Brown for (PetscInt d = 0; d < coord_dim; d++) { 11905f34f2dcSJed Brown for (PetscInt n = 0; n < num_local_nodes; n++) { 11915f34f2dcSJed Brown PetscInt gn = node_l2g[n]; 11925f34f2dcSJed Brown if (gn < nStart || nEnd <= gn) continue; 11935f34f2dcSJed Brown x[gn - nStart] = X[n * coord_dim + d]; 11945f34f2dcSJed Brown } 11955f34f2dcSJed Brown // CGNS nodes use 1-based indexing 11965f34f2dcSJed Brown cgsize_t start = nStart + 1, end = nEnd; 11974c155f28SJames Wright PetscCallCGNSWriteData(cgp_coord_write_data(cgv->file_num, base, zone, coord_ids[d], &start, &end, x), dm, viewer); 11985f34f2dcSJed Brown } 11995f34f2dcSJed Brown PetscCall(PetscFree(x)); 12005f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(coord, &X)); 12015f34f2dcSJed Brown } 12025f34f2dcSJed Brown 120359191b1eSJames Wright e_owned = cEnd - cStart; 120444698e90SJames Wright if (e_owned > 0) { 120544698e90SJames Wright DMPolytopeType cell_type; 120644698e90SJames Wright 120744698e90SJames Wright PetscCall(DMPlexGetCellType(dm, cStart, &cell_type)); 12085f34f2dcSJed Brown for (PetscInt i = cStart, c = 0; i < cEnd; i++) { 12095f34f2dcSJed Brown PetscInt closure_dof, *closure_indices, elem_size; 121044698e90SJames Wright 12115f34f2dcSJed Brown PetscCall(DMPlexGetClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 12125f34f2dcSJed Brown elem_size = closure_dof / coord_dim; 121344698e90SJames Wright if (!conn) PetscCall(PetscMalloc1(e_owned * elem_size, &conn)); 12145f34f2dcSJed Brown PetscCall(DMPlexCGNSGetPermutation_Internal(cell_type, closure_dof / coord_dim, &element_type, &perm)); 1215ad540459SPierre Jolivet for (PetscInt j = 0; j < elem_size; j++) conn[c++] = node_l2g[closure_indices[perm[j] * coord_dim] / coord_dim] + 1; 12165f34f2dcSJed Brown PetscCall(DMPlexRestoreClosureIndices(cdm, cdm->localSection, cdm->localSection, i, PETSC_FALSE, &closure_dof, &closure_indices, NULL, NULL)); 12175f34f2dcSJed Brown } 121844698e90SJames Wright } 121944698e90SJames Wright 122059191b1eSJames Wright { // Get global element_type (for ranks that do not have owned elements) 122159191b1eSJames Wright PetscInt local_element_type, global_element_type; 122259191b1eSJames Wright 122359191b1eSJames Wright local_element_type = e_owned > 0 ? element_type : -1; 1224462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&local_element_type, &global_element_type, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 122559191b1eSJames 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"); 1226ca328833SJames Wright element_type = (CGNS_ENUMT(ElementType_t))global_element_type; 122759191b1eSJames Wright } 1228462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&e_owned, &e_global, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 12291fb9530eSJeongu 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); 12305f34f2dcSJed Brown e_start = 0; 12311fb9530eSJeongu Kim PetscCallMPI(MPI_Exscan(&e_owned, &e_start, 1, MPIU_CGSIZE, MPI_SUM, PetscObjectComm((PetscObject)dm))); 12324c155f28SJames Wright PetscCallCGNSWrite(cgp_section_write(cgv->file_num, base, zone, "Elem", element_type, 1, e_global, 0, §ion), dm, viewer); 12334c155f28SJames Wright PetscCallCGNSWriteData(cgp_elements_write_data(cgv->file_num, base, zone, section, e_start + 1, e_start + e_owned, conn), dm, viewer); 12345f34f2dcSJed Brown PetscCall(PetscFree(conn)); 12355f34f2dcSJed Brown 12365f34f2dcSJed Brown cgv->base = base; 12375f34f2dcSJed Brown cgv->zone = zone; 12385f34f2dcSJed Brown cgv->node_l2g = node_l2g; 12395f34f2dcSJed Brown cgv->num_local_nodes = num_local_nodes; 12405f34f2dcSJed Brown cgv->nStart = nStart; 12415f34f2dcSJed Brown cgv->nEnd = nEnd; 12429bb8f83bSJed Brown cgv->eStart = e_start; 12439bb8f83bSJed Brown cgv->eEnd = e_start + e_owned; 12445f34f2dcSJed Brown if (1) { 12455f34f2dcSJed Brown PetscMPIInt rank; 12465f34f2dcSJed Brown int *efield; 12475f34f2dcSJed Brown int sol, field; 12485f34f2dcSJed Brown DMLabel label; 12495f34f2dcSJed Brown PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank)); 12505f34f2dcSJed Brown PetscCall(PetscMalloc1(e_owned, &efield)); 12515f34f2dcSJed Brown for (PetscInt i = 0; i < e_owned; i++) efield[i] = rank; 12524c155f28SJames Wright PetscCallCGNSWrite(cg_sol_write(cgv->file_num, base, zone, "CellInfo", CGNS_ENUMV(CellCenter), &sol), dm, viewer); 12534c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "Rank", &field), dm, viewer); 12545f34f2dcSJed Brown cgsize_t start = e_start + 1, end = e_start + e_owned; 12554c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 12565f34f2dcSJed Brown PetscCall(DMGetLabel(dm, "Cell Sets", &label)); 12575f34f2dcSJed Brown if (label) { 12585f34f2dcSJed Brown for (PetscInt c = cStart; c < cEnd; c++) { 12595f34f2dcSJed Brown PetscInt value; 12605f34f2dcSJed Brown PetscCall(DMLabelGetValue(label, c, &value)); 12615f34f2dcSJed Brown efield[c - cStart] = value; 12625f34f2dcSJed Brown } 12634c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, base, zone, sol, CGNS_ENUMV(Integer), "CellSet", &field), dm, viewer); 12644c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, base, zone, sol, field, &start, &end, efield), dm, viewer); 12655f34f2dcSJed Brown } 12665f34f2dcSJed Brown PetscCall(PetscFree(efield)); 12675f34f2dcSJed Brown } 12685f34f2dcSJed Brown } 1269b85bf5edSJed Brown PetscCall(DMDestroy(&colloc_dm)); 12703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12715f34f2dcSJed Brown } 12725f34f2dcSJed Brown 1273d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_Plex_Local_CGNS(Vec V, PetscViewer viewer) 1274d71ae5a4SJacob Faibussowitsch { 12755f34f2dcSJed Brown PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 12765f34f2dcSJed Brown DM dm; 12775f34f2dcSJed Brown PetscSection section; 127844698e90SJames Wright PetscInt time_step, num_fields, pStart, pEnd, fvGhostStart; 12795f34f2dcSJed Brown PetscReal time, *time_slot; 12809812b6beSJed Brown size_t *step_slot; 12815f34f2dcSJed Brown const PetscScalar *v; 12825f34f2dcSJed Brown char solution_name[PETSC_MAX_PATH_LEN]; 12835f34f2dcSJed Brown int sol; 12845f34f2dcSJed Brown 12855f34f2dcSJed Brown PetscFunctionBegin; 12865f34f2dcSJed Brown PetscCall(VecGetDM(V, &dm)); 1287ec2db9e4SJames Wright PetscCall(DMGetLocalSection(dm, §ion)); 1288ec2db9e4SJames Wright PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 128944698e90SJames Wright PetscCall(DMPlexGetCellTypeStratum(dm, DM_POLYTOPE_FV_GHOST, &fvGhostStart, NULL)); 129044698e90SJames Wright if (fvGhostStart >= 0) pEnd = fvGhostStart; 1291ec2db9e4SJames Wright 12925f34f2dcSJed Brown if (!cgv->node_l2g) PetscCall(DMView(dm, viewer)); 1293ec2db9e4SJames Wright if (!cgv->grid_loc) { // Determine if writing to cell-centers or to nodes 1294ec2db9e4SJames Wright PetscInt cStart, cEnd; 1295ec2db9e4SJames Wright PetscInt local_grid_loc, global_grid_loc; 1296ec2db9e4SJames Wright 1297ec2db9e4SJames Wright PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 129844698e90SJames Wright if (fvGhostStart >= 0) cEnd = fvGhostStart; 1299ec2db9e4SJames Wright if (cgv->num_local_nodes == 0) local_grid_loc = -1; 1300ec2db9e4SJames Wright else if (cStart == pStart && cEnd == pEnd) local_grid_loc = CGNS_ENUMV(CellCenter); 1301ec2db9e4SJames Wright else local_grid_loc = CGNS_ENUMV(Vertex); 1302ec2db9e4SJames Wright 1303462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&local_grid_loc, &global_grid_loc, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)viewer))); 1304ec2db9e4SJames Wright if (local_grid_loc != -1) 1305ec2db9e4SJames 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); 1306ec2db9e4SJames 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); 1307ca328833SJames Wright cgv->grid_loc = (CGNS_ENUMT(GridLocation_t))global_grid_loc; 1308ec2db9e4SJames Wright } 1309ec2db9e4SJames Wright if (!cgv->nodal_field) { 1310ec2db9e4SJames Wright switch (cgv->grid_loc) { 1311ec2db9e4SJames Wright case CGNS_ENUMV(Vertex): { 1312ec2db9e4SJames Wright PetscCall(PetscMalloc1(cgv->nEnd - cgv->nStart, &cgv->nodal_field)); 1313ec2db9e4SJames Wright } break; 1314ec2db9e4SJames Wright case CGNS_ENUMV(CellCenter): { 1315ec2db9e4SJames Wright PetscCall(PetscMalloc1(cgv->eEnd - cgv->eStart, &cgv->nodal_field)); 1316ec2db9e4SJames Wright } break; 1317ec2db9e4SJames Wright default: 1318ec2db9e4SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 1319ec2db9e4SJames Wright } 1320ec2db9e4SJames Wright } 13215f34f2dcSJed Brown if (!cgv->output_times) PetscCall(PetscSegBufferCreate(sizeof(PetscReal), 20, &cgv->output_times)); 13229812b6beSJed Brown if (!cgv->output_steps) PetscCall(PetscSegBufferCreate(sizeof(size_t), 20, &cgv->output_steps)); 13235f34f2dcSJed Brown 13245f34f2dcSJed Brown PetscCall(DMGetOutputSequenceNumber(dm, &time_step, &time)); 13259371c9d4SSatish Balay if (time_step < 0) { 13269371c9d4SSatish Balay time_step = 0; 13279371c9d4SSatish Balay time = 0.; 13289371c9d4SSatish Balay } 1329dfa0de3dSJames Wright // Avoid "Duplicate child name found" error by not writing an already-written solution. 1330dfa0de3dSJames Wright // This usually occurs when a solution is written and then diverges on the very next timestep. 1331dfa0de3dSJames Wright if (time_step == cgv->previous_output_step) PetscFunctionReturn(PETSC_SUCCESS); 1332dfa0de3dSJames Wright 13335f34f2dcSJed Brown PetscCall(PetscSegBufferGet(cgv->output_times, 1, &time_slot)); 13345f34f2dcSJed Brown *time_slot = time; 13359812b6beSJed Brown PetscCall(PetscSegBufferGet(cgv->output_steps, 1, &step_slot)); 1336dfa0de3dSJames Wright *step_slot = cgv->previous_output_step = time_step; 13375f34f2dcSJed Brown PetscCall(PetscSNPrintf(solution_name, sizeof solution_name, "FlowSolution%" PetscInt_FMT, time_step)); 13384c155f28SJames Wright PetscCallCGNSWrite(cg_sol_write(cgv->file_num, cgv->base, cgv->zone, solution_name, cgv->grid_loc, &sol), V, viewer); 13395f34f2dcSJed Brown PetscCall(VecGetArrayRead(V, &v)); 134000a86070SJed Brown PetscCall(PetscSectionGetNumFields(section, &num_fields)); 134100a86070SJed Brown for (PetscInt field = 0; field < num_fields; field++) { 134200a86070SJed Brown PetscInt ncomp; 13439bb8f83bSJed Brown const char *field_name; 13449bb8f83bSJed Brown PetscCall(PetscSectionGetFieldName(section, field, &field_name)); 134500a86070SJed Brown PetscCall(PetscSectionGetFieldComponents(section, field, &ncomp)); 13465f34f2dcSJed Brown for (PetscInt comp = 0; comp < ncomp; comp++) { 13475f34f2dcSJed Brown int cgfield; 13485f34f2dcSJed Brown const char *comp_name; 13499bb8f83bSJed Brown char cgns_field_name[32]; // CGNS max field name is 32 13505f34f2dcSJed Brown CGNS_ENUMT(DataType_t) datatype; 13515f34f2dcSJed Brown PetscCall(PetscSectionGetComponentName(section, field, comp, &comp_name)); 135244698e90SJames 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)); 13539bb8f83bSJed Brown else if (field_name[0] == '\0') PetscCall(PetscStrncpy(cgns_field_name, comp_name, sizeof cgns_field_name)); 13549bb8f83bSJed Brown else PetscCall(PetscSNPrintf(cgns_field_name, sizeof cgns_field_name, "%s.%s", field_name, comp_name)); 13555f34f2dcSJed Brown PetscCall(PetscCGNSDataType(PETSC_SCALAR, &datatype)); 13564c155f28SJames Wright PetscCallCGNSWrite(cgp_field_write(cgv->file_num, cgv->base, cgv->zone, sol, datatype, cgns_field_name, &cgfield), V, viewer); 135700a86070SJed Brown for (PetscInt p = pStart, n = 0; p < pEnd; p++) { 135800a86070SJed Brown PetscInt off, dof; 135900a86070SJed Brown PetscCall(PetscSectionGetFieldDof(section, p, field, &dof)); 136000a86070SJed Brown if (dof == 0) continue; 136100a86070SJed Brown PetscCall(PetscSectionGetFieldOffset(section, p, field, &off)); 136200a86070SJed Brown for (PetscInt c = comp; c < dof; c += ncomp, n++) { 1363ec2db9e4SJames Wright switch (cgv->grid_loc) { 13649bb8f83bSJed Brown case CGNS_ENUMV(Vertex): { 13655f34f2dcSJed Brown PetscInt gn = cgv->node_l2g[n]; 13665f34f2dcSJed Brown if (gn < cgv->nStart || cgv->nEnd <= gn) continue; 136700a86070SJed Brown cgv->nodal_field[gn - cgv->nStart] = v[off + c]; 13689bb8f83bSJed Brown } break; 13699bb8f83bSJed Brown case CGNS_ENUMV(CellCenter): { 13709bb8f83bSJed Brown cgv->nodal_field[n] = v[off + c]; 13719bb8f83bSJed Brown } break; 13729bb8f83bSJed Brown default: 13739bb8f83bSJed Brown SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only pack for Vertex and CellCenter grid locations"); 13749bb8f83bSJed Brown } 137500a86070SJed Brown } 13765f34f2dcSJed Brown } 13775f34f2dcSJed Brown // CGNS nodes use 1-based indexing 1378ec2db9e4SJames Wright cgsize_t start, end; 1379ec2db9e4SJames Wright switch (cgv->grid_loc) { 1380ec2db9e4SJames Wright case CGNS_ENUMV(Vertex): { 1381ec2db9e4SJames Wright start = cgv->nStart + 1; 1382ec2db9e4SJames Wright end = cgv->nEnd; 1383ec2db9e4SJames Wright } break; 1384ec2db9e4SJames Wright case CGNS_ENUMV(CellCenter): { 13859bb8f83bSJed Brown start = cgv->eStart + 1; 13869bb8f83bSJed Brown end = cgv->eEnd; 1387ec2db9e4SJames Wright } break; 1388ec2db9e4SJames Wright default: 1389ec2db9e4SJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Can only write for Vertex and CellCenter grid locations"); 13909bb8f83bSJed Brown } 13914c155f28SJames Wright PetscCallCGNSWriteData(cgp_field_write_data(cgv->file_num, cgv->base, cgv->zone, sol, cgfield, &start, &end, cgv->nodal_field), V, viewer); 13925f34f2dcSJed Brown } 139300a86070SJed Brown } 13945f34f2dcSJed Brown PetscCall(VecRestoreArrayRead(V, &v)); 139525760affSJed Brown PetscCall(PetscViewerCGNSCheckBatch_Internal(viewer)); 13963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13975f34f2dcSJed Brown } 1398472b9844SJames Wright 1399472b9844SJames Wright PetscErrorCode VecLoad_Plex_CGNS_Internal(Vec V, PetscViewer viewer) 1400472b9844SJames Wright { 1401472b9844SJames Wright MPI_Comm comm; 1402472b9844SJames Wright char buffer[CGIO_MAX_NAME_LENGTH + 1]; 1403472b9844SJames Wright PetscViewer_CGNS *cgv = (PetscViewer_CGNS *)viewer->data; 1404472b9844SJames Wright int cgid = cgv->file_num; 1405472b9844SJames Wright PetscBool use_parallel_viewer = PETSC_FALSE; 1406472b9844SJames Wright int z = 1; // Only support one zone 1407472b9844SJames Wright int B = 1; // Only support one base 1408472b9844SJames Wright int numComp; 1409472b9844SJames Wright PetscInt V_numComps, mystartv, myendv, myownedv; 1410472b9844SJames Wright 1411b30057acSJames Wright PetscFunctionBegin; 1412472b9844SJames Wright PetscCall(PetscObjectGetComm((PetscObject)V, &comm)); 1413472b9844SJames Wright 1414472b9844SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_plex_cgns_parallel", &use_parallel_viewer, NULL)); 1415472b9844SJames 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"); 1416472b9844SJames Wright 1417472b9844SJames Wright { // Get CGNS node ownership information 1418472b9844SJames Wright int nbases, nzones; 1419472b9844SJames Wright PetscInt NVertices; 1420472b9844SJames Wright PetscLayout vtx_map; 1421472b9844SJames Wright cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 1422472b9844SJames Wright 14234c155f28SJames Wright PetscCallCGNSRead(cg_nbases(cgid, &nbases), V, viewer); 1424472b9844SJames Wright PetscCheck(nbases <= 1, PETSC_COMM_SELF, PETSC_ERR_LIB, "CGNS file must have a single base, not %d", nbases); 14254c155f28SJames Wright PetscCallCGNSRead(cg_nzones(cgid, B, &nzones), V, viewer); 1426472b9844SJames Wright PetscCheck(nzones == 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "limited to one zone %d", (int)nzones); 1427472b9844SJames Wright 14284c155f28SJames Wright PetscCallCGNSRead(cg_zone_read(cgid, B, z, buffer, sizes), V, viewer); 1429472b9844SJames Wright NVertices = sizes[0]; 1430472b9844SJames Wright 1431472b9844SJames Wright PetscCall(PetscLayoutCreateFromSizes(comm, PETSC_DECIDE, NVertices, 1, &vtx_map)); 1432472b9844SJames Wright PetscCall(PetscLayoutGetRange(vtx_map, &mystartv, &myendv)); 1433472b9844SJames Wright PetscCall(PetscLayoutGetLocalSize(vtx_map, &myownedv)); 1434472b9844SJames Wright PetscCall(PetscLayoutDestroy(&vtx_map)); 1435472b9844SJames Wright } 1436472b9844SJames Wright 1437472b9844SJames Wright { // -- Read data from file into Vec 1438472b9844SJames Wright PetscScalar *fields = NULL; 1439472b9844SJames Wright PetscSF sfNatural; 1440472b9844SJames Wright 1441472b9844SJames Wright { // Check compatibility between sfNatural and the data source and sink 1442472b9844SJames Wright DM dm; 1443472b9844SJames Wright PetscInt nleaves, nroots, V_local_size; 1444472b9844SJames Wright 1445472b9844SJames Wright PetscCall(VecGetDM(V, &dm)); 1446472b9844SJames Wright PetscCall(DMGetNaturalSF(dm, &sfNatural)); 1447472b9844SJames Wright PetscCheck(sfNatural, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "DM of Vec must have sfNatural"); 1448472b9844SJames Wright PetscCall(PetscSFGetGraph(sfNatural, &nroots, &nleaves, NULL, NULL)); 1449472b9844SJames Wright PetscCall(VecGetLocalSize(V, &V_local_size)); 1450472b9844SJames 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); 1451*d0bbc01eSJames Wright if (nroots == 0) { 1452*d0bbc01eSJames Wright PetscCheck(V_local_size == nroots, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Local Vec size (% " PetscInt_FMT ") must be zero if number of roots in sfNatural is zero", V_local_size); 1453*d0bbc01eSJames Wright V_numComps = 0; 1454*d0bbc01eSJames Wright } else { 1455472b9844SJames 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); 1456472b9844SJames Wright V_numComps = V_local_size / nroots; 1457472b9844SJames Wright } 1458*d0bbc01eSJames Wright } 1459472b9844SJames Wright 1460472b9844SJames Wright { // Read data into component-major ordering 1461472b9844SJames Wright int isol, numSols; 1462472b9844SJames Wright CGNS_ENUMT(DataType_t) datatype; 1463472b9844SJames Wright double *fields_CGNS; 1464472b9844SJames Wright 14654c155f28SJames Wright PetscCallCGNSRead(cg_nsols(cgid, B, z, &numSols), V, viewer); 1466472b9844SJames Wright PetscCall(PetscViewerCGNSGetSolutionFileIndex_Internal(viewer, &isol)); 14674c155f28SJames Wright PetscCallCGNSRead(cg_nfields(cgid, B, z, isol, &numComp), V, viewer); 1468*d0bbc01eSJames Wright PetscCheck(V_numComps == numComp || V_numComps == 0, 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); 1469472b9844SJames Wright 1470472b9844SJames Wright cgsize_t range_min[3] = {mystartv + 1, 1, 1}; 1471472b9844SJames Wright cgsize_t range_max[3] = {myendv, 1, 1}; 1472472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * numComp, &fields_CGNS)); 1473472b9844SJames Wright PetscCall(PetscMalloc1(myownedv * numComp, &fields)); 1474472b9844SJames Wright for (int d = 0; d < numComp; ++d) { 14754c155f28SJames Wright PetscCallCGNSRead(cg_field_info(cgid, B, z, isol, (d + 1), &datatype, buffer), V, viewer); 1476472b9844SJames Wright PetscCheck(datatype == CGNS_ENUMV(RealDouble), PETSC_COMM_SELF, PETSC_ERR_ARG_NOTSAMETYPE, "Field %s in file is not of type double", buffer); 14774c155f28SJames Wright PetscCallCGNSReadData(cgp_field_read_data(cgid, B, z, isol, (d + 1), range_min, range_max, &fields_CGNS[d * myownedv]), V, viewer); 1478472b9844SJames Wright } 1479472b9844SJames Wright for (int d = 0; d < numComp; ++d) { 1480472b9844SJames Wright for (PetscInt v = 0; v < myownedv; ++v) fields[v * numComp + d] = fields_CGNS[d * myownedv + v]; 1481472b9844SJames Wright } 1482472b9844SJames Wright PetscCall(PetscFree(fields_CGNS)); 1483472b9844SJames Wright } 1484472b9844SJames Wright 1485472b9844SJames Wright { // Reduce fields into Vec array 1486472b9844SJames Wright PetscScalar *V_array; 1487472b9844SJames Wright MPI_Datatype fieldtype; 1488472b9844SJames Wright 1489472b9844SJames Wright PetscCall(VecGetArrayWrite(V, &V_array)); 1490472b9844SJames Wright PetscCallMPI(MPI_Type_contiguous(numComp, MPIU_SCALAR, &fieldtype)); 1491472b9844SJames Wright PetscCallMPI(MPI_Type_commit(&fieldtype)); 1492472b9844SJames Wright PetscCall(PetscSFReduceBegin(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1493472b9844SJames Wright PetscCall(PetscSFReduceEnd(sfNatural, fieldtype, fields, V_array, MPI_REPLACE)); 1494472b9844SJames Wright PetscCallMPI(MPI_Type_free(&fieldtype)); 1495472b9844SJames Wright PetscCall(VecRestoreArrayWrite(V, &V_array)); 1496472b9844SJames Wright } 1497472b9844SJames Wright PetscCall(PetscFree(fields)); 1498472b9844SJames Wright } 1499472b9844SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1500472b9844SJames Wright } 1501