1e1b39ce3SMatthew G. Knepley #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3e1b39ce3SMatthew G. Knepley 41b7ec9c1SMatthew G. Knepley #undef I /* Very old CGNS stupidly uses I as a variable, which fails when using complex. Curse you idiot package managers */ 5e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 6e1b39ce3SMatthew G. Knepley #include <cgnslib.h> 7e1b39ce3SMatthew G. Knepley #include <cgns_io.h> 8e1b39ce3SMatthew G. Knepley #endif 9c4cbe8f1SToby Isaac #if !defined(CGNS_ENUMT) 10c4cbe8f1SToby Isaac #define CGNS_ENUMT(a) a 11c4cbe8f1SToby Isaac #endif 12c4cbe8f1SToby Isaac #if !defined(CGNS_ENUMV) 13c4cbe8f1SToby Isaac #define CGNS_ENUMV(a) a 14c4cbe8f1SToby Isaac #endif 15e1b39ce3SMatthew G. Knepley 16c261946bSLisandro Dalcin #define CHKERRCGNS(ierr) \ 17c261946bSLisandro Dalcin do { \ 18c261946bSLisandro Dalcin int _cgns_ier = (ierr); \ 19*28b400f6SJacob Faibussowitsch PetscCheck(!_cgns_ier,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS error %d %s",_cgns_ier,cg_get_error()); \ 20c261946bSLisandro Dalcin } while (0) 21c261946bSLisandro Dalcin 2244cd5272SMichael Lange /*@C 2344cd5272SMichael Lange DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file. 2444cd5272SMichael Lange 25d083f849SBarry Smith Collective 2644cd5272SMichael Lange 2744cd5272SMichael Lange Input Parameters: 2844cd5272SMichael Lange + comm - The MPI communicator 2944cd5272SMichael Lange . filename - The name of the CGNS file 3044cd5272SMichael Lange - interpolate - Create faces and edges in the mesh 3144cd5272SMichael Lange 3244cd5272SMichael Lange Output Parameter: 3344cd5272SMichael Lange . dm - The DM object representing the mesh 3444cd5272SMichael Lange 3544cd5272SMichael Lange Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html 3644cd5272SMichael Lange 3744cd5272SMichael Lange Level: beginner 3844cd5272SMichael Lange 3944cd5272SMichael Lange .seealso: DMPlexCreate(), DMPlexCreateCGNS(), DMPlexCreateExodus() 4044cd5272SMichael Lange @*/ 4144cd5272SMichael Lange PetscErrorCode DMPlexCreateCGNSFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 4244cd5272SMichael Lange { 4344cd5272SMichael Lange PetscMPIInt rank; 4444cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS) 4544cd5272SMichael Lange int cgid = -1; 4644cd5272SMichael Lange #endif 4744cd5272SMichael Lange 4844cd5272SMichael Lange PetscFunctionBegin; 4944cd5272SMichael Lange PetscValidCharPointer(filename, 2); 505f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 5144cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS) 52dd400576SPatrick Sanan if (rank == 0) { 535f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_open(filename, CG_MODE_READ, &cgid)); 542c71b3e2SJacob Faibussowitsch PetscCheckFalse(cgid <= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 5544cd5272SMichael Lange } 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexCreateCGNS(comm, cgid, interpolate, dm)); 575f80ce2aSJacob Faibussowitsch if (rank == 0) CHKERRCGNS(cg_close(cgid)); 58b458e8f1SJose E. Roman PetscFunctionReturn(0); 5944cd5272SMichael Lange #else 6044cd5272SMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir"); 6144cd5272SMichael Lange #endif 6244cd5272SMichael Lange } 6344cd5272SMichael Lange 64e1b39ce3SMatthew G. Knepley /*@ 6544cd5272SMichael Lange DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID. 66e1b39ce3SMatthew G. Knepley 67d083f849SBarry Smith Collective 68e1b39ce3SMatthew G. Knepley 69e1b39ce3SMatthew G. Knepley Input Parameters: 70e1b39ce3SMatthew G. Knepley + comm - The MPI communicator 71e1b39ce3SMatthew G. Knepley . cgid - The CG id associated with a file and obtained using cg_open 72e1b39ce3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 73e1b39ce3SMatthew G. Knepley 74e1b39ce3SMatthew G. Knepley Output Parameter: 75e1b39ce3SMatthew G. Knepley . dm - The DM object representing the mesh 76e1b39ce3SMatthew G. Knepley 77e1b39ce3SMatthew G. Knepley Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html 78e1b39ce3SMatthew G. Knepley 79e1b39ce3SMatthew G. Knepley Level: beginner 80e1b39ce3SMatthew G. Knepley 81e1b39ce3SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexCreateExodus() 82e1b39ce3SMatthew G. Knepley @*/ 83e1b39ce3SMatthew G. Knepley PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 84e1b39ce3SMatthew G. Knepley { 85e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 86e1b39ce3SMatthew G. Knepley PetscMPIInt num_proc, rank; 87c261946bSLisandro Dalcin DM cdm; 884e6cc61cSLisandro Dalcin DMLabel label; 89e1b39ce3SMatthew G. Knepley PetscSection coordSection; 90e1b39ce3SMatthew G. Knepley Vec coordinates; 91e1b39ce3SMatthew G. Knepley PetscScalar *coords; 92c261946bSLisandro Dalcin PetscInt *cellStart, *vertStart, v; 934e6cc61cSLisandro Dalcin PetscInt labelIdRange[2], labelId; 94e1b39ce3SMatthew G. Knepley PetscErrorCode ierr; 95e1b39ce3SMatthew G. Knepley /* Read from file */ 96e1b39ce3SMatthew G. Knepley char basename[CGIO_MAX_NAME_LENGTH+1]; 97e1b39ce3SMatthew G. Knepley char buffer[CGIO_MAX_NAME_LENGTH+1]; 98c261946bSLisandro Dalcin int dim = 0, physDim = 0, coordDim =0, numVertices = 0, numCells = 0; 99e1b39ce3SMatthew G. Knepley int nzones = 0; 100e1b39ce3SMatthew G. Knepley #endif 101e1b39ce3SMatthew G. Knepley 102e1b39ce3SMatthew G. Knepley PetscFunctionBegin; 103e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 1045f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 1055f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &num_proc)); 1065f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreate(comm, dm)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetType(*dm, DMPLEX)); 108c261946bSLisandro Dalcin 109a5b23f4aSJose E. Roman /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 110dd400576SPatrick Sanan if (rank == 0) { 111e1b39ce3SMatthew G. Knepley int nbases, z; 112e1b39ce3SMatthew G. Knepley 1135f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_nbases(cgid, &nbases)); 1142c71b3e2SJacob Faibussowitsch PetscCheckFalse(nbases > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d",nbases); 1155f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_base_read(cgid, 1, basename, &dim, &physDim)); 1165f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_nzones(cgid, 1, &nzones)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart)); 118e1b39ce3SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 119e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 120e1b39ce3SMatthew G. Knepley 1215f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_zone_read(cgid, 1, z, buffer, sizes)); 122e1b39ce3SMatthew G. Knepley numVertices += sizes[0]; 123e1b39ce3SMatthew G. Knepley numCells += sizes[1]; 124f8716870SMatthew G. Knepley cellStart[z] += sizes[1] + cellStart[z-1]; 125f8716870SMatthew G. Knepley vertStart[z] += sizes[0] + vertStart[z-1]; 126f8716870SMatthew G. Knepley } 127f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 128f8716870SMatthew G. Knepley vertStart[z] += numCells; 129e1b39ce3SMatthew G. Knepley } 130c261946bSLisandro Dalcin coordDim = dim; 131e1b39ce3SMatthew G. Knepley } 1325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm)); 1335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 1345f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm)); 1355f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&nzones, 1, MPI_INT, 0, comm)); 136c261946bSLisandro Dalcin 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) *dm, basename)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetDimension(*dm, dim)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(*dm, "celltype")); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetChart(*dm, 0, numCells+numVertices)); 141e1b39ce3SMatthew G. Knepley 142e1b39ce3SMatthew G. Knepley /* Read zone information */ 143dd400576SPatrick Sanan if (rank == 0) { 144c261946bSLisandro Dalcin int z, c, c_loc; 145e1b39ce3SMatthew G. Knepley 146e1b39ce3SMatthew G. Knepley /* Read the cell set connectivity table and build mesh topology 147e1b39ce3SMatthew G. Knepley CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 148e1b39ce3SMatthew G. Knepley /* First set sizes */ 149e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 150c4cbe8f1SToby Isaac CGNS_ENUMT(ZoneType_t) zonetype; 151e1b39ce3SMatthew G. Knepley int nsections; 152c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 153e1b39ce3SMatthew G. Knepley cgsize_t start, end; 154e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 155e1b39ce3SMatthew G. Knepley PetscInt numCorners; 156c261946bSLisandro Dalcin DMPolytopeType ctype; 157e1b39ce3SMatthew G. Knepley 1585f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_zone_type(cgid, 1, z, &zonetype)); 1592c71b3e2SJacob Faibussowitsch PetscCheckFalse(zonetype == CGNS_ENUMV(Structured),PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS"); 1605f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_nsections(cgid, 1, z, &nsections)); 1612c71b3e2SJacob Faibussowitsch PetscCheckFalse(nsections > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d",nsections); 1625f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag)); 163e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 164c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 165e1b39ce3SMatthew G. Knepley cgsize_t elementDataSize, *elements; 166e1b39ce3SMatthew G. Knepley PetscInt off; 167e1b39ce3SMatthew G. Knepley 1685f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(elementDataSize, &elements)); 1705f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL)); 17159dfc837SStefan Lemvig Glimberg for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 1724c9a562dSMatthew G. Knepley switch (elements[off]) { 173c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 174c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 175c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 176c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 177c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 178c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 17998921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]); 180e1b39ce3SMatthew G. Knepley } 1815f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetConeSize(*dm, c, numCorners)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCellType(*dm, c, ctype)); 183e1b39ce3SMatthew G. Knepley off += numCorners+1; 184e1b39ce3SMatthew G. Knepley } 1855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(elements)); 186e1b39ce3SMatthew G. Knepley } else { 187e1b39ce3SMatthew G. Knepley switch (cellType) { 188c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 189c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 190c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 191c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 192c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 193c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 19498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 195e1b39ce3SMatthew G. Knepley } 19659dfc837SStefan Lemvig Glimberg for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 1975f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetConeSize(*dm, c, numCorners)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCellType(*dm, c, ctype)); 199e1b39ce3SMatthew G. Knepley } 200e1b39ce3SMatthew G. Knepley } 201e1b39ce3SMatthew G. Knepley } 202c261946bSLisandro Dalcin for (v = numCells; v < numCells+numVertices; ++v) { 2035f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 204c261946bSLisandro Dalcin } 205c261946bSLisandro Dalcin } 206c261946bSLisandro Dalcin 2075f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(*dm)); 208c261946bSLisandro Dalcin 2095f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(*dm, "zone")); 210dd400576SPatrick Sanan if (rank == 0) { 211c261946bSLisandro Dalcin int z, c, c_loc, v_loc; 212c261946bSLisandro Dalcin 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, "zone", &label)); 214e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 215c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 216c261946bSLisandro Dalcin cgsize_t elementDataSize, *elements, start, end; 217e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 218e1b39ce3SMatthew G. Knepley PetscInt *cone, numc, numCorners, maxCorners = 27; 219e1b39ce3SMatthew G. Knepley 2205f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag)); 221e1b39ce3SMatthew G. Knepley numc = end - start; 222e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 2235f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(elementDataSize,&elements,maxCorners,&cone)); 2255f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL)); 226c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 227eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 22859dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 2294c9a562dSMatthew G. Knepley switch (elements[v]) { 230c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 231c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 232c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 233c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 234c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 235c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 23698921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]); 237e1b39ce3SMatthew G. Knepley } 2384c9a562dSMatthew G. Knepley ++v; 239e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 240e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 241e1b39ce3SMatthew G. Knepley } 2425f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexReorderCell(*dm, c, cone)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCone(*dm, c, cone)); 2445f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(label, c, z)); 245e1b39ce3SMatthew G. Knepley } 246e1b39ce3SMatthew G. Knepley } else { 247e1b39ce3SMatthew G. Knepley switch (cellType) { 248c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 249c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 250c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 251c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 252c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 253c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 25498921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 255e1b39ce3SMatthew G. Knepley } 256eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 25759dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 258e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 259e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 260e1b39ce3SMatthew G. Knepley } 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexReorderCell(*dm, c, cone)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSetCone(*dm, c, cone)); 2635f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(label, c, z)); 264e1b39ce3SMatthew G. Knepley } 265e1b39ce3SMatthew G. Knepley } 2665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(elements,cone)); 267e1b39ce3SMatthew G. Knepley } 268e1b39ce3SMatthew G. Knepley } 269c261946bSLisandro Dalcin 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexSymmetrize(*dm)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexStratify(*dm)); 272e1b39ce3SMatthew G. Knepley if (interpolate) { 2735fd9971aSMatthew G. Knepley DM idm; 274e1b39ce3SMatthew G. Knepley 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMPlexInterpolate(*dm, &idm)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(dm)); 277e1b39ce3SMatthew G. Knepley *dm = idm; 278e1b39ce3SMatthew G. Knepley } 279e1b39ce3SMatthew G. Knepley 280e1b39ce3SMatthew G. Knepley /* Read coordinates */ 2815f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinateDim(*dm, coordDim)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(*dm, &cdm)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalSection(cdm, &coordSection)); 2845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetNumFields(coordSection, 1)); 2855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldComponents(coordSection, 0, coordDim)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 287e1b39ce3SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) { 2885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetDof(coordSection, v, dim)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetFieldDof(coordSection, v, 0, coordDim)); 290e1b39ce3SMatthew G. Knepley } 2915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSectionSetUp(coordSection)); 292c261946bSLisandro Dalcin 2935f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLocalVector(cdm, &coordinates)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(coordinates, &coords)); 295dd400576SPatrick Sanan if (rank == 0) { 296e1b39ce3SMatthew G. Knepley PetscInt off = 0; 297e1b39ce3SMatthew G. Knepley float *x[3]; 29859dfc837SStefan Lemvig Glimberg int z, d; 299e1b39ce3SMatthew G. Knepley 3005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2])); 30159dfc837SStefan Lemvig Glimberg for (z = 1; z <= nzones; ++z) { 302c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 303e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 304e1b39ce3SMatthew G. Knepley cgsize_t range_min[3] = {1, 1, 1}; 305e1b39ce3SMatthew G. Knepley cgsize_t range_max[3] = {1, 1, 1}; 306e1b39ce3SMatthew G. Knepley int ngrids, ncoords; 307e1b39ce3SMatthew G. Knepley 3085f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_zone_read(cgid, 1, z, buffer, sizes)); 309e1b39ce3SMatthew G. Knepley range_max[0] = sizes[0]; 3105f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_ngrids(cgid, 1, z, &ngrids)); 3112c71b3e2SJacob Faibussowitsch PetscCheckFalse(ngrids > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d",ngrids); 3125f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_ncoords(cgid, 1, z, &ncoords)); 3132c71b3e2SJacob Faibussowitsch PetscCheckFalse(ncoords != coordDim,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d",ncoords); 314c261946bSLisandro Dalcin for (d = 0; d < coordDim; ++d) { 3155f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer)); 3165f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d])); 317e1b39ce3SMatthew G. Knepley } 318c261946bSLisandro Dalcin if (coordDim >= 1) { 319c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v]; 320e1b39ce3SMatthew G. Knepley } 321c261946bSLisandro Dalcin if (coordDim >= 2) { 322c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v]; 323e1b39ce3SMatthew G. Knepley } 324c261946bSLisandro Dalcin if (coordDim >= 3) { 325c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v]; 326e1b39ce3SMatthew G. Knepley } 327e1b39ce3SMatthew G. Knepley off += sizes[0]; 328e1b39ce3SMatthew G. Knepley } 3295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(x[0],x[1],x[2])); 330e1b39ce3SMatthew G. Knepley } 3315f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(coordinates, &coords)); 332c261946bSLisandro Dalcin 3335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 3345f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetBlockSize(coordinates, coordDim)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetCoordinatesLocal(*dm, coordinates)); 3365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&coordinates)); 337c261946bSLisandro Dalcin 338f8716870SMatthew G. Knepley /* Read boundary conditions */ 3395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumLabels(*dm, &labelIdRange[0])); 340dd400576SPatrick Sanan if (rank == 0) { 341c4cbe8f1SToby Isaac CGNS_ENUMT(BCType_t) bctype; 342c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 343c4cbe8f1SToby Isaac CGNS_ENUMT(PointSetType_t) pointtype; 344f8716870SMatthew G. Knepley cgsize_t *points; 345f8716870SMatthew G. Knepley PetscReal *normals; 346f8716870SMatthew G. Knepley int normal[3]; 347f8716870SMatthew G. Knepley char *bcname = buffer; 348f8716870SMatthew G. Knepley cgsize_t npoints, nnormals; 349f8716870SMatthew G. Knepley int z, nbc, bc, c, ndatasets; 350f8716870SMatthew G. Knepley 351f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 3525f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_nbocos(cgid, 1, z, &nbc)); 353f8716870SMatthew G. Knepley for (bc = 1; bc <= nbc; ++bc) { 3545f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets)); 3555f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateLabel(*dm, bcname)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabel(*dm, bcname, &label)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(npoints, &points, nnormals, &normals)); 3585f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_boco_read(cgid, 1, z, bc, points, (void *) normals)); 359c4cbe8f1SToby Isaac if (pointtype == CGNS_ENUMV(ElementRange)) { 360f8716870SMatthew G. Knepley /* Range of cells: assuming half-open interval since the documentation sucks */ 361f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 3625f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(label, c - cellStart[z-1], 1)); 363f8716870SMatthew G. Knepley } 364c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(ElementList)) { 365f8716870SMatthew G. Knepley /* List of cells */ 366f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 3675f80ce2aSJacob Faibussowitsch CHKERRQ(DMLabelSetValue(label, points[c] - cellStart[z-1], 1)); 368f8716870SMatthew G. Knepley } 369c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointRange)) { 370c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 371f8716870SMatthew G. Knepley 372f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 3735f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 3745f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_gridlocation_read(&gridloc)); 375f8716870SMatthew G. Knepley /* Range of points: assuming half-open interval since the documentation sucks */ 376f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 3775f80ce2aSJacob Faibussowitsch if (gridloc == CGNS_ENUMV(Vertex)) CHKERRQ(DMLabelSetValue(label, c - vertStart[z-1], 1)); 3785f80ce2aSJacob Faibussowitsch else CHKERRQ(DMLabelSetValue(label, c - cellStart[z-1], 1)); 379f8716870SMatthew G. Knepley } 380c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointList)) { 381c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 382f8716870SMatthew G. Knepley 383f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 3845f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end")); 3855f80ce2aSJacob Faibussowitsch CHKERRCGNS(cg_gridlocation_read(&gridloc)); 386f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 3875f80ce2aSJacob Faibussowitsch if (gridloc == CGNS_ENUMV(Vertex)) CHKERRQ(DMLabelSetValue(label, points[c] - vertStart[z-1], 1)); 3885f80ce2aSJacob Faibussowitsch else CHKERRQ(DMLabelSetValue(label, points[c] - cellStart[z-1], 1)); 389f8716870SMatthew G. Knepley } 39098921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype); 3915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(points, normals)); 392f8716870SMatthew G. Knepley } 393f8716870SMatthew G. Knepley } 3945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(cellStart, vertStart)); 395f8716870SMatthew G. Knepley } 3965f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetNumLabels(*dm, &labelIdRange[1])); 3975f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm)); 3984e6cc61cSLisandro Dalcin 3994e6cc61cSLisandro Dalcin /* Create BC labels at all processes */ 4004e6cc61cSLisandro Dalcin for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 4014e6cc61cSLisandro Dalcin char *labelName = buffer; 4024e6cc61cSLisandro Dalcin size_t len = sizeof(buffer); 4034e6cc61cSLisandro Dalcin const char *locName; 4044e6cc61cSLisandro Dalcin 405dd400576SPatrick Sanan if (rank == 0) { 4065f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLabelByNum(*dm, labelId, &label)); 4075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetName((PetscObject)label, &locName)); 4085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(labelName, locName, len)); 4094e6cc61cSLisandro Dalcin } 4105f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm)); 4115f80ce2aSJacob Faibussowitsch CHKERRMPI(DMCreateLabel(*dm, labelName)); 4124e6cc61cSLisandro Dalcin } 413b458e8f1SJose E. Roman PetscFunctionReturn(0); 414e1b39ce3SMatthew G. Knepley #else 415e1b39ce3SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir"); 416e1b39ce3SMatthew G. Knepley #endif 417e1b39ce3SMatthew G. Knepley } 418