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 16*c261946bSLisandro Dalcin #define CHKERRCGNS(ierr) \ 17*c261946bSLisandro Dalcin do { \ 18*c261946bSLisandro Dalcin int _cgns_ier = (ierr); \ 19*c261946bSLisandro Dalcin if (PetscUnlikely(_cgns_ier)) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS error %d %s",_cgns_ier,cg_get_error()); \ 20*c261946bSLisandro Dalcin } while (0) 21*c261946bSLisandro 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 PetscErrorCode ierr; 4544cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS) 4644cd5272SMichael Lange int cgid = -1; 4744cd5272SMichael Lange #endif 4844cd5272SMichael Lange 4944cd5272SMichael Lange PetscFunctionBegin; 5044cd5272SMichael Lange PetscValidCharPointer(filename, 2); 51ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 5244cd5272SMichael Lange #if defined(PETSC_HAVE_CGNS) 5344cd5272SMichael Lange if (!rank) { 54*c261946bSLisandro Dalcin ierr = cg_open(filename, CG_MODE_READ, &cgid);CHKERRCGNS(ierr); 5544cd5272SMichael Lange if (cgid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "cg_open(\"%s\",...) did not return a valid file ID", filename); 5644cd5272SMichael Lange } 5744cd5272SMichael Lange ierr = DMPlexCreateCGNS(comm, cgid, interpolate, dm);CHKERRQ(ierr); 58*c261946bSLisandro Dalcin if (!rank) {ierr = cg_close(cgid);CHKERRCGNS(ierr);} 59b458e8f1SJose E. Roman PetscFunctionReturn(0); 6044cd5272SMichael Lange #else 6144cd5272SMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "Loading meshes requires CGNS support. Reconfigure using --with-cgns-dir"); 6244cd5272SMichael Lange #endif 6344cd5272SMichael Lange } 6444cd5272SMichael Lange 65e1b39ce3SMatthew G. Knepley /*@ 6644cd5272SMichael Lange DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file ID. 67e1b39ce3SMatthew G. Knepley 68d083f849SBarry Smith Collective 69e1b39ce3SMatthew G. Knepley 70e1b39ce3SMatthew G. Knepley Input Parameters: 71e1b39ce3SMatthew G. Knepley + comm - The MPI communicator 72e1b39ce3SMatthew G. Knepley . cgid - The CG id associated with a file and obtained using cg_open 73e1b39ce3SMatthew G. Knepley - interpolate - Create faces and edges in the mesh 74e1b39ce3SMatthew G. Knepley 75e1b39ce3SMatthew G. Knepley Output Parameter: 76e1b39ce3SMatthew G. Knepley . dm - The DM object representing the mesh 77e1b39ce3SMatthew G. Knepley 78e1b39ce3SMatthew G. Knepley Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html 79e1b39ce3SMatthew G. Knepley 80e1b39ce3SMatthew G. Knepley Level: beginner 81e1b39ce3SMatthew G. Knepley 82e1b39ce3SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexCreateExodus() 83e1b39ce3SMatthew G. Knepley @*/ 84e1b39ce3SMatthew G. Knepley PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 85e1b39ce3SMatthew G. Knepley { 86e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 87e1b39ce3SMatthew G. Knepley PetscMPIInt num_proc, rank; 88*c261946bSLisandro Dalcin DM cdm; 89e1b39ce3SMatthew G. Knepley PetscSection coordSection; 90e1b39ce3SMatthew G. Knepley Vec coordinates; 91e1b39ce3SMatthew G. Knepley PetscScalar *coords; 92*c261946bSLisandro Dalcin PetscInt *cellStart, *vertStart, v; 93e1b39ce3SMatthew G. Knepley PetscErrorCode ierr; 94e1b39ce3SMatthew G. Knepley /* Read from file */ 95e1b39ce3SMatthew G. Knepley char basename[CGIO_MAX_NAME_LENGTH+1]; 96e1b39ce3SMatthew G. Knepley char buffer[CGIO_MAX_NAME_LENGTH+1]; 97*c261946bSLisandro Dalcin int dim = 0, physDim = 0, coordDim =0, numVertices = 0, numCells = 0; 98e1b39ce3SMatthew G. Knepley int nzones = 0; 99e1b39ce3SMatthew G. Knepley #endif 100e1b39ce3SMatthew G. Knepley 101e1b39ce3SMatthew G. Knepley PetscFunctionBegin; 102e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 103ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 104ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &num_proc);CHKERRMPI(ierr); 105e1b39ce3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 106e1b39ce3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 107*c261946bSLisandro Dalcin 108e1b39ce3SMatthew G. Knepley /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */ 109e1b39ce3SMatthew G. Knepley if (!rank) { 110e1b39ce3SMatthew G. Knepley int nbases, z; 111e1b39ce3SMatthew G. Knepley 112*c261946bSLisandro Dalcin ierr = cg_nbases(cgid, &nbases);CHKERRCGNS(ierr); 113e1b39ce3SMatthew G. Knepley if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases); 114*c261946bSLisandro Dalcin ierr = cg_base_read(cgid, 1, basename, &dim, &physDim);CHKERRCGNS(ierr); 115*c261946bSLisandro Dalcin ierr = cg_nzones(cgid, 1, &nzones);CHKERRCGNS(ierr); 116f8716870SMatthew G. Knepley ierr = PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart);CHKERRQ(ierr); 117e1b39ce3SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 118e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 119e1b39ce3SMatthew G. Knepley 120*c261946bSLisandro Dalcin ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRCGNS(ierr); 121e1b39ce3SMatthew G. Knepley numVertices += sizes[0]; 122e1b39ce3SMatthew G. Knepley numCells += sizes[1]; 123f8716870SMatthew G. Knepley cellStart[z] += sizes[1] + cellStart[z-1]; 124f8716870SMatthew G. Knepley vertStart[z] += sizes[0] + vertStart[z-1]; 125f8716870SMatthew G. Knepley } 126f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 127f8716870SMatthew G. Knepley vertStart[z] += numCells; 128e1b39ce3SMatthew G. Knepley } 129*c261946bSLisandro Dalcin coordDim = dim; 130e1b39ce3SMatthew G. Knepley } 131ffc4695bSBarry Smith ierr = MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);CHKERRMPI(ierr); 132ffc4695bSBarry Smith ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 133*c261946bSLisandro Dalcin ierr = MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 134ffc4695bSBarry Smith ierr = MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 135*c261946bSLisandro Dalcin 136e1b39ce3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, basename);CHKERRQ(ierr); 137c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 138*c261946bSLisandro Dalcin ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr); 139e1b39ce3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 140e1b39ce3SMatthew G. Knepley 141e1b39ce3SMatthew G. Knepley /* Read zone information */ 142e1b39ce3SMatthew G. Knepley if (!rank) { 143*c261946bSLisandro Dalcin int z, c, c_loc; 144e1b39ce3SMatthew G. Knepley 145e1b39ce3SMatthew G. Knepley /* Read the cell set connectivity table and build mesh topology 146e1b39ce3SMatthew G. Knepley CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 147e1b39ce3SMatthew G. Knepley /* First set sizes */ 148e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 149c4cbe8f1SToby Isaac CGNS_ENUMT(ZoneType_t) zonetype; 150e1b39ce3SMatthew G. Knepley int nsections; 151c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 152e1b39ce3SMatthew G. Knepley cgsize_t start, end; 153e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 154e1b39ce3SMatthew G. Knepley PetscInt numCorners; 155*c261946bSLisandro Dalcin DMPolytopeType ctype; 156e1b39ce3SMatthew G. Knepley 157*c261946bSLisandro Dalcin ierr = cg_zone_type(cgid, 1, z, &zonetype);CHKERRCGNS(ierr); 158c4cbe8f1SToby Isaac if (zonetype == CGNS_ENUMV(Structured)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS"); 159*c261946bSLisandro Dalcin ierr = cg_nsections(cgid, 1, z, &nsections);CHKERRCGNS(ierr); 160e1b39ce3SMatthew G. Knepley if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections); 161*c261946bSLisandro Dalcin ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRCGNS(ierr); 162e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 163c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 164e1b39ce3SMatthew G. Knepley cgsize_t elementDataSize, *elements; 165e1b39ce3SMatthew G. Knepley PetscInt off; 166e1b39ce3SMatthew G. Knepley 167*c261946bSLisandro Dalcin ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRCGNS(ierr); 168785e854fSJed Brown ierr = PetscMalloc1(elementDataSize, &elements);CHKERRQ(ierr); 169*c261946bSLisandro Dalcin ierr = cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL);CHKERRCGNS(ierr); 17059dfc837SStefan Lemvig Glimberg for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 1714c9a562dSMatthew G. Knepley switch (elements[off]) { 172*c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 173*c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 174*c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 175*c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 176*c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 177*c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 1784c9a562dSMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]); 179e1b39ce3SMatthew G. Knepley } 180e1b39ce3SMatthew G. Knepley ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 181*c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, c, ctype);CHKERRQ(ierr); 182e1b39ce3SMatthew G. Knepley off += numCorners+1; 183e1b39ce3SMatthew G. Knepley } 184e1b39ce3SMatthew G. Knepley ierr = PetscFree(elements);CHKERRQ(ierr); 185e1b39ce3SMatthew G. Knepley } else { 186e1b39ce3SMatthew G. Knepley switch (cellType) { 187*c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 188*c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 189*c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 190*c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 191*c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 192*c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 193e1b39ce3SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 194e1b39ce3SMatthew G. Knepley } 19559dfc837SStefan Lemvig Glimberg for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 196e1b39ce3SMatthew G. Knepley ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 197*c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, c, ctype);CHKERRQ(ierr); 198e1b39ce3SMatthew G. Knepley } 199e1b39ce3SMatthew G. Knepley } 200e1b39ce3SMatthew G. Knepley } 201*c261946bSLisandro Dalcin for (v = numCells; v < numCells+numVertices; ++v) { 202*c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr); 203*c261946bSLisandro Dalcin } 204*c261946bSLisandro Dalcin } 205*c261946bSLisandro Dalcin 206e1b39ce3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 207*c261946bSLisandro Dalcin 208*c261946bSLisandro Dalcin if (!rank) { 209*c261946bSLisandro Dalcin int z, c, c_loc, v_loc; 210*c261946bSLisandro Dalcin 211e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 212c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 213*c261946bSLisandro Dalcin cgsize_t elementDataSize, *elements, start, end; 214e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 215e1b39ce3SMatthew G. Knepley PetscInt *cone, numc, numCorners, maxCorners = 27; 216e1b39ce3SMatthew G. Knepley 217*c261946bSLisandro Dalcin ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRCGNS(ierr); 218e1b39ce3SMatthew G. Knepley numc = end - start; 219e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 220*c261946bSLisandro Dalcin ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRCGNS(ierr); 221dcca6d9dSJed Brown ierr = PetscMalloc2(elementDataSize,&elements,maxCorners,&cone);CHKERRQ(ierr); 222*c261946bSLisandro Dalcin ierr = cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL);CHKERRCGNS(ierr); 223c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 224eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 22559dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 2264c9a562dSMatthew G. Knepley switch (elements[v]) { 227*c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 228c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 229c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 230c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 231*c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 232c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 2334c9a562dSMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]); 234e1b39ce3SMatthew G. Knepley } 2354c9a562dSMatthew G. Knepley ++v; 236e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 237e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 238e1b39ce3SMatthew G. Knepley } 239*c261946bSLisandro Dalcin ierr = DMPlexReorderCell(*dm, c, cone);CHKERRQ(ierr); 240e1b39ce3SMatthew G. Knepley ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 241c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr); 242e1b39ce3SMatthew G. Knepley } 243e1b39ce3SMatthew G. Knepley } else { 244e1b39ce3SMatthew G. Knepley switch (cellType) { 245*c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 246c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 247c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 248c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 249*c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 250c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 251e1b39ce3SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 252e1b39ce3SMatthew G. Knepley } 253eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 25459dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 255e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 256e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 257e1b39ce3SMatthew G. Knepley } 258*c261946bSLisandro Dalcin ierr = DMPlexReorderCell(*dm, c, cone);CHKERRQ(ierr); 259e1b39ce3SMatthew G. Knepley ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 260c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr); 261e1b39ce3SMatthew G. Knepley } 262e1b39ce3SMatthew G. Knepley } 263e1b39ce3SMatthew G. Knepley ierr = PetscFree2(elements,cone);CHKERRQ(ierr); 264e1b39ce3SMatthew G. Knepley } 265e1b39ce3SMatthew G. Knepley } 266*c261946bSLisandro Dalcin 267e1b39ce3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 268e1b39ce3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 269e1b39ce3SMatthew G. Knepley if (interpolate) { 2705fd9971aSMatthew G. Knepley DM idm; 271e1b39ce3SMatthew G. Knepley 272e1b39ce3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 273e1b39ce3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 274e1b39ce3SMatthew G. Knepley *dm = idm; 275e1b39ce3SMatthew G. Knepley } 276e1b39ce3SMatthew G. Knepley 277e1b39ce3SMatthew G. Knepley /* Read coordinates */ 278*c261946bSLisandro Dalcin ierr = DMSetCoordinateDim(*dm, coordDim);CHKERRQ(ierr); 279*c261946bSLisandro Dalcin ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 280*c261946bSLisandro Dalcin ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 281e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 282*c261946bSLisandro Dalcin ierr = PetscSectionSetFieldComponents(coordSection, 0, coordDim);CHKERRQ(ierr); 283e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 284e1b39ce3SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) { 285e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 286*c261946bSLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, v, 0, coordDim);CHKERRQ(ierr); 287e1b39ce3SMatthew G. Knepley } 288e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 289*c261946bSLisandro Dalcin 290*c261946bSLisandro Dalcin ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 291e1b39ce3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 292e1b39ce3SMatthew G. Knepley if (!rank) { 293e1b39ce3SMatthew G. Knepley PetscInt off = 0; 294e1b39ce3SMatthew G. Knepley float *x[3]; 29559dfc837SStefan Lemvig Glimberg int z, d; 296e1b39ce3SMatthew G. Knepley 297dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]);CHKERRQ(ierr); 29859dfc837SStefan Lemvig Glimberg for (z = 1; z <= nzones; ++z) { 299c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 300e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 301e1b39ce3SMatthew G. Knepley cgsize_t range_min[3] = {1, 1, 1}; 302e1b39ce3SMatthew G. Knepley cgsize_t range_max[3] = {1, 1, 1}; 303e1b39ce3SMatthew G. Knepley int ngrids, ncoords; 304e1b39ce3SMatthew G. Knepley 305*c261946bSLisandro Dalcin ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRCGNS(ierr); 306e1b39ce3SMatthew G. Knepley range_max[0] = sizes[0]; 307*c261946bSLisandro Dalcin ierr = cg_ngrids(cgid, 1, z, &ngrids);CHKERRCGNS(ierr); 308e1b39ce3SMatthew G. Knepley if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids); 309*c261946bSLisandro Dalcin ierr = cg_ncoords(cgid, 1, z, &ncoords);CHKERRCGNS(ierr); 310*c261946bSLisandro Dalcin if (ncoords != coordDim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d\n",ncoords); 311*c261946bSLisandro Dalcin for (d = 0; d < coordDim; ++d) { 312*c261946bSLisandro Dalcin ierr = cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer);CHKERRCGNS(ierr); 313*c261946bSLisandro Dalcin ierr = cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]);CHKERRCGNS(ierr); 314e1b39ce3SMatthew G. Knepley } 315*c261946bSLisandro Dalcin if (coordDim >= 1) { 316*c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v]; 317e1b39ce3SMatthew G. Knepley } 318*c261946bSLisandro Dalcin if (coordDim >= 2) { 319*c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v]; 320e1b39ce3SMatthew G. Knepley } 321*c261946bSLisandro Dalcin if (coordDim >= 3) { 322*c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v]; 323e1b39ce3SMatthew G. Knepley } 324e1b39ce3SMatthew G. Knepley off += sizes[0]; 325e1b39ce3SMatthew G. Knepley } 326e1b39ce3SMatthew G. Knepley ierr = PetscFree3(x[0],x[1],x[2]);CHKERRQ(ierr); 327e1b39ce3SMatthew G. Knepley } 328e1b39ce3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 329*c261946bSLisandro Dalcin 330*c261946bSLisandro Dalcin ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 331*c261946bSLisandro Dalcin ierr = VecSetBlockSize(coordinates, coordDim);CHKERRQ(ierr); 332e1b39ce3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 333e1b39ce3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 334*c261946bSLisandro Dalcin 335f8716870SMatthew G. Knepley /* Read boundary conditions */ 336f8716870SMatthew G. Knepley if (!rank) { 337f8716870SMatthew G. Knepley DMLabel label; 338c4cbe8f1SToby Isaac CGNS_ENUMT(BCType_t) bctype; 339c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 340c4cbe8f1SToby Isaac CGNS_ENUMT(PointSetType_t) pointtype; 341f8716870SMatthew G. Knepley cgsize_t *points; 342f8716870SMatthew G. Knepley PetscReal *normals; 343f8716870SMatthew G. Knepley int normal[3]; 344f8716870SMatthew G. Knepley char *bcname = buffer; 345f8716870SMatthew G. Knepley cgsize_t npoints, nnormals; 346f8716870SMatthew G. Knepley int z, nbc, bc, c, ndatasets; 347f8716870SMatthew G. Knepley 348f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 349*c261946bSLisandro Dalcin ierr = cg_nbocos(cgid, 1, z, &nbc);CHKERRCGNS(ierr); 350f8716870SMatthew G. Knepley for (bc = 1; bc <= nbc; ++bc) { 351*c261946bSLisandro Dalcin ierr = cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets);CHKERRCGNS(ierr); 352c58f1c22SToby Isaac ierr = DMCreateLabel(*dm, bcname);CHKERRQ(ierr); 353c58f1c22SToby Isaac ierr = DMGetLabel(*dm, bcname, &label);CHKERRQ(ierr); 354f8716870SMatthew G. Knepley ierr = PetscMalloc2(npoints, &points, nnormals, &normals);CHKERRQ(ierr); 355*c261946bSLisandro Dalcin ierr = cg_boco_read(cgid, 1, z, bc, points, (void *) normals);CHKERRCGNS(ierr); 356c4cbe8f1SToby Isaac if (pointtype == CGNS_ENUMV(ElementRange)) { 357f8716870SMatthew G. Knepley /* Range of cells: assuming half-open interval since the documentation sucks */ 358f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 359f8716870SMatthew G. Knepley ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr); 360f8716870SMatthew G. Knepley } 361c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(ElementList)) { 362f8716870SMatthew G. Knepley /* List of cells */ 363f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 364f8716870SMatthew G. Knepley ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr); 365f8716870SMatthew G. Knepley } 366c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointRange)) { 367c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 368f8716870SMatthew G. Knepley 369f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 370*c261946bSLisandro Dalcin ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");CHKERRCGNS(ierr); 371*c261946bSLisandro Dalcin ierr = cg_gridlocation_read(&gridloc);CHKERRCGNS(ierr); 372f8716870SMatthew G. Knepley /* Range of points: assuming half-open interval since the documentation sucks */ 373f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 374c4cbe8f1SToby Isaac if (gridloc == CGNS_ENUMV(Vertex)) {ierr = DMLabelSetValue(label, c - vertStart[z-1], 1);CHKERRQ(ierr);} 375f8716870SMatthew G. Knepley else {ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr);} 376f8716870SMatthew G. Knepley } 377c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointList)) { 378c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 379f8716870SMatthew G. Knepley 380f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 381*c261946bSLisandro Dalcin ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");CHKERRCGNS(ierr); 382*c261946bSLisandro Dalcin ierr = cg_gridlocation_read(&gridloc);CHKERRCGNS(ierr); 383f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 384c4cbe8f1SToby Isaac if (gridloc == CGNS_ENUMV(Vertex)) {ierr = DMLabelSetValue(label, points[c] - vertStart[z-1], 1);CHKERRQ(ierr);} 385f8716870SMatthew G. Knepley else {ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr);} 386f8716870SMatthew G. Knepley } 387f8716870SMatthew G. Knepley } else SETERRQ1(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype); 388f8716870SMatthew G. Knepley ierr = PetscFree2(points, normals);CHKERRQ(ierr); 389f8716870SMatthew G. Knepley } 390f8716870SMatthew G. Knepley } 391f8716870SMatthew G. Knepley ierr = PetscFree2(cellStart, vertStart);CHKERRQ(ierr); 392f8716870SMatthew G. Knepley } 393b458e8f1SJose E. Roman PetscFunctionReturn(0); 394e1b39ce3SMatthew G. Knepley #else 395e1b39ce3SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir"); 396e1b39ce3SMatthew G. Knepley #endif 397e1b39ce3SMatthew G. Knepley } 398