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*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(_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 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) 53dd400576SPatrick Sanan if (rank == 0) { 54c261946bSLisandro Dalcin ierr = cg_open(filename, CG_MODE_READ, &cgid);CHKERRCGNS(ierr); 55*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(cgid <= 0,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); 58dd400576SPatrick Sanan if (rank == 0) {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; 88c261946bSLisandro Dalcin DM cdm; 894e6cc61cSLisandro Dalcin DMLabel label; 90e1b39ce3SMatthew G. Knepley PetscSection coordSection; 91e1b39ce3SMatthew G. Knepley Vec coordinates; 92e1b39ce3SMatthew G. Knepley PetscScalar *coords; 93c261946bSLisandro Dalcin PetscInt *cellStart, *vertStart, v; 944e6cc61cSLisandro Dalcin PetscInt labelIdRange[2], labelId; 95e1b39ce3SMatthew G. Knepley PetscErrorCode ierr; 96e1b39ce3SMatthew G. Knepley /* Read from file */ 97e1b39ce3SMatthew G. Knepley char basename[CGIO_MAX_NAME_LENGTH+1]; 98e1b39ce3SMatthew G. Knepley char buffer[CGIO_MAX_NAME_LENGTH+1]; 99c261946bSLisandro Dalcin int dim = 0, physDim = 0, coordDim =0, numVertices = 0, numCells = 0; 100e1b39ce3SMatthew G. Knepley int nzones = 0; 101e1b39ce3SMatthew G. Knepley #endif 102e1b39ce3SMatthew G. Knepley 103e1b39ce3SMatthew G. Knepley PetscFunctionBegin; 104e1b39ce3SMatthew G. Knepley #if defined(PETSC_HAVE_CGNS) 105ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr); 106ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &num_proc);CHKERRMPI(ierr); 107e1b39ce3SMatthew G. Knepley ierr = DMCreate(comm, dm);CHKERRQ(ierr); 108e1b39ce3SMatthew G. Knepley ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 109c261946bSLisandro Dalcin 110a5b23f4aSJose E. Roman /* Open CGNS II file and read basic information on rank 0, then broadcast to all processors */ 111dd400576SPatrick Sanan if (rank == 0) { 112e1b39ce3SMatthew G. Knepley int nbases, z; 113e1b39ce3SMatthew G. Knepley 114c261946bSLisandro Dalcin ierr = cg_nbases(cgid, &nbases);CHKERRCGNS(ierr); 115*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nbases > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d",nbases); 116c261946bSLisandro Dalcin ierr = cg_base_read(cgid, 1, basename, &dim, &physDim);CHKERRCGNS(ierr); 117c261946bSLisandro Dalcin ierr = cg_nzones(cgid, 1, &nzones);CHKERRCGNS(ierr); 118f8716870SMatthew G. Knepley ierr = PetscCalloc2(nzones+1, &cellStart, nzones+1, &vertStart);CHKERRQ(ierr); 119e1b39ce3SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 120e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 121e1b39ce3SMatthew G. Knepley 122c261946bSLisandro Dalcin ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRCGNS(ierr); 123e1b39ce3SMatthew G. Knepley numVertices += sizes[0]; 124e1b39ce3SMatthew G. Knepley numCells += sizes[1]; 125f8716870SMatthew G. Knepley cellStart[z] += sizes[1] + cellStart[z-1]; 126f8716870SMatthew G. Knepley vertStart[z] += sizes[0] + vertStart[z-1]; 127f8716870SMatthew G. Knepley } 128f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 129f8716870SMatthew G. Knepley vertStart[z] += numCells; 130e1b39ce3SMatthew G. Knepley } 131c261946bSLisandro Dalcin coordDim = dim; 132e1b39ce3SMatthew G. Knepley } 133ffc4695bSBarry Smith ierr = MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);CHKERRMPI(ierr); 134ffc4695bSBarry Smith ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 135c261946bSLisandro Dalcin ierr = MPI_Bcast(&coordDim, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 136ffc4695bSBarry Smith ierr = MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);CHKERRMPI(ierr); 137c261946bSLisandro Dalcin 138e1b39ce3SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dm, basename);CHKERRQ(ierr); 139c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 140c261946bSLisandro Dalcin ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr); 141e1b39ce3SMatthew G. Knepley ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 142e1b39ce3SMatthew G. Knepley 143e1b39ce3SMatthew G. Knepley /* Read zone information */ 144dd400576SPatrick Sanan if (rank == 0) { 145c261946bSLisandro Dalcin int z, c, c_loc; 146e1b39ce3SMatthew G. Knepley 147e1b39ce3SMatthew G. Knepley /* Read the cell set connectivity table and build mesh topology 148e1b39ce3SMatthew G. Knepley CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 149e1b39ce3SMatthew G. Knepley /* First set sizes */ 150e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 151c4cbe8f1SToby Isaac CGNS_ENUMT(ZoneType_t) zonetype; 152e1b39ce3SMatthew G. Knepley int nsections; 153c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 154e1b39ce3SMatthew G. Knepley cgsize_t start, end; 155e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 156e1b39ce3SMatthew G. Knepley PetscInt numCorners; 157c261946bSLisandro Dalcin DMPolytopeType ctype; 158e1b39ce3SMatthew G. Knepley 159c261946bSLisandro Dalcin ierr = cg_zone_type(cgid, 1, z, &zonetype);CHKERRCGNS(ierr); 160*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(zonetype == CGNS_ENUMV(Structured),PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS"); 161c261946bSLisandro Dalcin ierr = cg_nsections(cgid, 1, z, &nsections);CHKERRCGNS(ierr); 162*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(nsections > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d",nsections); 163c261946bSLisandro Dalcin ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRCGNS(ierr); 164e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 165c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 166e1b39ce3SMatthew G. Knepley cgsize_t elementDataSize, *elements; 167e1b39ce3SMatthew G. Knepley PetscInt off; 168e1b39ce3SMatthew G. Knepley 169c261946bSLisandro Dalcin ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRCGNS(ierr); 170785e854fSJed Brown ierr = PetscMalloc1(elementDataSize, &elements);CHKERRQ(ierr); 171c261946bSLisandro Dalcin ierr = cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL);CHKERRCGNS(ierr); 17259dfc837SStefan Lemvig Glimberg for (c_loc = start, off = 0; c_loc <= end; ++c_loc, ++c) { 1734c9a562dSMatthew G. Knepley switch (elements[off]) { 174c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 175c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 176c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 177c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 178c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 179c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 18098921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[off]); 181e1b39ce3SMatthew G. Knepley } 182e1b39ce3SMatthew G. Knepley ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 183c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, c, ctype);CHKERRQ(ierr); 184e1b39ce3SMatthew G. Knepley off += numCorners+1; 185e1b39ce3SMatthew G. Knepley } 186e1b39ce3SMatthew G. Knepley ierr = PetscFree(elements);CHKERRQ(ierr); 187e1b39ce3SMatthew G. Knepley } else { 188e1b39ce3SMatthew G. Knepley switch (cellType) { 189c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; ctype = DM_POLYTOPE_SEGMENT; break; 190c261946bSLisandro Dalcin case CGNS_ENUMV(TRI_3): numCorners = 3; ctype = DM_POLYTOPE_TRIANGLE; break; 191c261946bSLisandro Dalcin case CGNS_ENUMV(QUAD_4): numCorners = 4; ctype = DM_POLYTOPE_QUADRILATERAL; break; 192c261946bSLisandro Dalcin case CGNS_ENUMV(TETRA_4): numCorners = 4; ctype = DM_POLYTOPE_TETRAHEDRON; break; 193c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; ctype = DM_POLYTOPE_TRI_PRISM; break; 194c261946bSLisandro Dalcin case CGNS_ENUMV(HEXA_8): numCorners = 8; ctype = DM_POLYTOPE_HEXAHEDRON; break; 19598921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 196e1b39ce3SMatthew G. Knepley } 19759dfc837SStefan Lemvig Glimberg for (c_loc = start; c_loc <= end; ++c_loc, ++c) { 198e1b39ce3SMatthew G. Knepley ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 199c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, c, ctype);CHKERRQ(ierr); 200e1b39ce3SMatthew G. Knepley } 201e1b39ce3SMatthew G. Knepley } 202e1b39ce3SMatthew G. Knepley } 203c261946bSLisandro Dalcin for (v = numCells; v < numCells+numVertices; ++v) { 204c261946bSLisandro Dalcin ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr); 205c261946bSLisandro Dalcin } 206c261946bSLisandro Dalcin } 207c261946bSLisandro Dalcin 208e1b39ce3SMatthew G. Knepley ierr = DMSetUp(*dm);CHKERRQ(ierr); 209c261946bSLisandro Dalcin 2104e6cc61cSLisandro Dalcin ierr = DMCreateLabel(*dm, "zone");CHKERRQ(ierr); 211dd400576SPatrick Sanan if (rank == 0) { 212c261946bSLisandro Dalcin int z, c, c_loc, v_loc; 213c261946bSLisandro Dalcin 2144e6cc61cSLisandro Dalcin ierr = DMGetLabel(*dm, "zone", &label);CHKERRQ(ierr); 215e1b39ce3SMatthew G. Knepley for (z = 1, c = 0; z <= nzones; ++z) { 216c4cbe8f1SToby Isaac CGNS_ENUMT(ElementType_t) cellType; 217c261946bSLisandro Dalcin cgsize_t elementDataSize, *elements, start, end; 218e1b39ce3SMatthew G. Knepley int nbndry, parentFlag; 219e1b39ce3SMatthew G. Knepley PetscInt *cone, numc, numCorners, maxCorners = 27; 220e1b39ce3SMatthew G. Knepley 221c261946bSLisandro Dalcin ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRCGNS(ierr); 222e1b39ce3SMatthew G. Knepley numc = end - start; 223e1b39ce3SMatthew G. Knepley /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 224c261946bSLisandro Dalcin ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRCGNS(ierr); 225dcca6d9dSJed Brown ierr = PetscMalloc2(elementDataSize,&elements,maxCorners,&cone);CHKERRQ(ierr); 226c261946bSLisandro Dalcin ierr = cg_poly_elements_read(cgid, 1, z, 1, elements, NULL, NULL);CHKERRCGNS(ierr); 227c4cbe8f1SToby Isaac if (cellType == CGNS_ENUMV(MIXED)) { 228eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 22959dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 2304c9a562dSMatthew G. Knepley switch (elements[v]) { 231c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 232c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 233c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 234c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 235c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 236c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 23798921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) elements[v]); 238e1b39ce3SMatthew G. Knepley } 2394c9a562dSMatthew G. Knepley ++v; 240e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 241e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 242e1b39ce3SMatthew G. Knepley } 243c261946bSLisandro Dalcin ierr = DMPlexReorderCell(*dm, c, cone);CHKERRQ(ierr); 244e1b39ce3SMatthew G. Knepley ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 2454e6cc61cSLisandro Dalcin ierr = DMLabelSetValue(label, c, z);CHKERRQ(ierr); 246e1b39ce3SMatthew G. Knepley } 247e1b39ce3SMatthew G. Knepley } else { 248e1b39ce3SMatthew G. Knepley switch (cellType) { 249c261946bSLisandro Dalcin case CGNS_ENUMV(BAR_2): numCorners = 2; break; 250c4cbe8f1SToby Isaac case CGNS_ENUMV(TRI_3): numCorners = 3; break; 251c4cbe8f1SToby Isaac case CGNS_ENUMV(QUAD_4): numCorners = 4; break; 252c4cbe8f1SToby Isaac case CGNS_ENUMV(TETRA_4): numCorners = 4; break; 253c261946bSLisandro Dalcin case CGNS_ENUMV(PENTA_6): numCorners = 6; break; 254c4cbe8f1SToby Isaac case CGNS_ENUMV(HEXA_8): numCorners = 8; break; 25598921bdaSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 256e1b39ce3SMatthew G. Knepley } 257eaf898f9SPatrick Sanan /* CGNS uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 25859dfc837SStefan Lemvig Glimberg for (c_loc = 0, v = 0; c_loc <= numc; ++c_loc, ++c) { 259e1b39ce3SMatthew G. Knepley for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 260e1b39ce3SMatthew G. Knepley cone[v_loc] = elements[v]+numCells-1; 261e1b39ce3SMatthew G. Knepley } 262c261946bSLisandro Dalcin ierr = DMPlexReorderCell(*dm, c, cone);CHKERRQ(ierr); 263e1b39ce3SMatthew G. Knepley ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 2644e6cc61cSLisandro Dalcin ierr = DMLabelSetValue(label, c, z);CHKERRQ(ierr); 265e1b39ce3SMatthew G. Knepley } 266e1b39ce3SMatthew G. Knepley } 267e1b39ce3SMatthew G. Knepley ierr = PetscFree2(elements,cone);CHKERRQ(ierr); 268e1b39ce3SMatthew G. Knepley } 269e1b39ce3SMatthew G. Knepley } 270c261946bSLisandro Dalcin 271e1b39ce3SMatthew G. Knepley ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 272e1b39ce3SMatthew G. Knepley ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 273e1b39ce3SMatthew G. Knepley if (interpolate) { 2745fd9971aSMatthew G. Knepley DM idm; 275e1b39ce3SMatthew G. Knepley 276e1b39ce3SMatthew G. Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 277e1b39ce3SMatthew G. Knepley ierr = DMDestroy(dm);CHKERRQ(ierr); 278e1b39ce3SMatthew G. Knepley *dm = idm; 279e1b39ce3SMatthew G. Knepley } 280e1b39ce3SMatthew G. Knepley 281e1b39ce3SMatthew G. Knepley /* Read coordinates */ 282c261946bSLisandro Dalcin ierr = DMSetCoordinateDim(*dm, coordDim);CHKERRQ(ierr); 283c261946bSLisandro Dalcin ierr = DMGetCoordinateDM(*dm, &cdm);CHKERRQ(ierr); 284c261946bSLisandro Dalcin ierr = DMGetLocalSection(cdm, &coordSection);CHKERRQ(ierr); 285e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 286c261946bSLisandro Dalcin ierr = PetscSectionSetFieldComponents(coordSection, 0, coordDim);CHKERRQ(ierr); 287e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 288e1b39ce3SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) { 289e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 290c261946bSLisandro Dalcin ierr = PetscSectionSetFieldDof(coordSection, v, 0, coordDim);CHKERRQ(ierr); 291e1b39ce3SMatthew G. Knepley } 292e1b39ce3SMatthew G. Knepley ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 293c261946bSLisandro Dalcin 294c261946bSLisandro Dalcin ierr = DMCreateLocalVector(cdm, &coordinates);CHKERRQ(ierr); 295e1b39ce3SMatthew G. Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 296dd400576SPatrick Sanan if (rank == 0) { 297e1b39ce3SMatthew G. Knepley PetscInt off = 0; 298e1b39ce3SMatthew G. Knepley float *x[3]; 29959dfc837SStefan Lemvig Glimberg int z, d; 300e1b39ce3SMatthew G. Knepley 301dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x[0],numVertices,&x[1],numVertices,&x[2]);CHKERRQ(ierr); 30259dfc837SStefan Lemvig Glimberg for (z = 1; z <= nzones; ++z) { 303c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 304e1b39ce3SMatthew G. Knepley cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 305e1b39ce3SMatthew G. Knepley cgsize_t range_min[3] = {1, 1, 1}; 306e1b39ce3SMatthew G. Knepley cgsize_t range_max[3] = {1, 1, 1}; 307e1b39ce3SMatthew G. Knepley int ngrids, ncoords; 308e1b39ce3SMatthew G. Knepley 309c261946bSLisandro Dalcin ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRCGNS(ierr); 310e1b39ce3SMatthew G. Knepley range_max[0] = sizes[0]; 311c261946bSLisandro Dalcin ierr = cg_ngrids(cgid, 1, z, &ngrids);CHKERRCGNS(ierr); 312*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(ngrids > 1,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d",ngrids); 313c261946bSLisandro Dalcin ierr = cg_ncoords(cgid, 1, z, &ncoords);CHKERRCGNS(ierr); 314*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(ncoords != coordDim,PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d",ncoords); 315c261946bSLisandro Dalcin for (d = 0; d < coordDim; ++d) { 316c261946bSLisandro Dalcin ierr = cg_coord_info(cgid, 1, z, 1+d, &datatype, buffer);CHKERRCGNS(ierr); 317c261946bSLisandro Dalcin ierr = cg_coord_read(cgid, 1, z, buffer, CGNS_ENUMV(RealSingle), range_min, range_max, x[d]);CHKERRCGNS(ierr); 318e1b39ce3SMatthew G. Knepley } 319c261946bSLisandro Dalcin if (coordDim >= 1) { 320c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+0] = x[0][v]; 321e1b39ce3SMatthew G. Knepley } 322c261946bSLisandro Dalcin if (coordDim >= 2) { 323c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+1] = x[1][v]; 324e1b39ce3SMatthew G. Knepley } 325c261946bSLisandro Dalcin if (coordDim >= 3) { 326c261946bSLisandro Dalcin for (v = 0; v < sizes[0]; ++v) coords[(v+off)*coordDim+2] = x[2][v]; 327e1b39ce3SMatthew G. Knepley } 328e1b39ce3SMatthew G. Knepley off += sizes[0]; 329e1b39ce3SMatthew G. Knepley } 330e1b39ce3SMatthew G. Knepley ierr = PetscFree3(x[0],x[1],x[2]);CHKERRQ(ierr); 331e1b39ce3SMatthew G. Knepley } 332e1b39ce3SMatthew G. Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 333c261946bSLisandro Dalcin 334c261946bSLisandro Dalcin ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 335c261946bSLisandro Dalcin ierr = VecSetBlockSize(coordinates, coordDim);CHKERRQ(ierr); 336e1b39ce3SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 337e1b39ce3SMatthew G. Knepley ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 338c261946bSLisandro Dalcin 339f8716870SMatthew G. Knepley /* Read boundary conditions */ 3404e6cc61cSLisandro Dalcin ierr = DMGetNumLabels(*dm, &labelIdRange[0]);CHKERRQ(ierr); 341dd400576SPatrick Sanan if (rank == 0) { 342c4cbe8f1SToby Isaac CGNS_ENUMT(BCType_t) bctype; 343c4cbe8f1SToby Isaac CGNS_ENUMT(DataType_t) datatype; 344c4cbe8f1SToby Isaac CGNS_ENUMT(PointSetType_t) pointtype; 345f8716870SMatthew G. Knepley cgsize_t *points; 346f8716870SMatthew G. Knepley PetscReal *normals; 347f8716870SMatthew G. Knepley int normal[3]; 348f8716870SMatthew G. Knepley char *bcname = buffer; 349f8716870SMatthew G. Knepley cgsize_t npoints, nnormals; 350f8716870SMatthew G. Knepley int z, nbc, bc, c, ndatasets; 351f8716870SMatthew G. Knepley 352f8716870SMatthew G. Knepley for (z = 1; z <= nzones; ++z) { 353c261946bSLisandro Dalcin ierr = cg_nbocos(cgid, 1, z, &nbc);CHKERRCGNS(ierr); 354f8716870SMatthew G. Knepley for (bc = 1; bc <= nbc; ++bc) { 355c261946bSLisandro Dalcin ierr = cg_boco_info(cgid, 1, z, bc, bcname, &bctype, &pointtype, &npoints, normal, &nnormals, &datatype, &ndatasets);CHKERRCGNS(ierr); 356c58f1c22SToby Isaac ierr = DMCreateLabel(*dm, bcname);CHKERRQ(ierr); 357c58f1c22SToby Isaac ierr = DMGetLabel(*dm, bcname, &label);CHKERRQ(ierr); 358f8716870SMatthew G. Knepley ierr = PetscMalloc2(npoints, &points, nnormals, &normals);CHKERRQ(ierr); 359c261946bSLisandro Dalcin ierr = cg_boco_read(cgid, 1, z, bc, points, (void *) normals);CHKERRCGNS(ierr); 360c4cbe8f1SToby Isaac if (pointtype == CGNS_ENUMV(ElementRange)) { 361f8716870SMatthew G. Knepley /* Range of cells: assuming half-open interval since the documentation sucks */ 362f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 363f8716870SMatthew G. Knepley ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr); 364f8716870SMatthew G. Knepley } 365c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(ElementList)) { 366f8716870SMatthew G. Knepley /* List of cells */ 367f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 368f8716870SMatthew G. Knepley ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr); 369f8716870SMatthew G. Knepley } 370c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointRange)) { 371c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 372f8716870SMatthew G. Knepley 373f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 374c261946bSLisandro Dalcin ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");CHKERRCGNS(ierr); 375c261946bSLisandro Dalcin ierr = cg_gridlocation_read(&gridloc);CHKERRCGNS(ierr); 376f8716870SMatthew G. Knepley /* Range of points: assuming half-open interval since the documentation sucks */ 377f8716870SMatthew G. Knepley for (c = points[0]; c < points[1]; ++c) { 378c4cbe8f1SToby Isaac if (gridloc == CGNS_ENUMV(Vertex)) {ierr = DMLabelSetValue(label, c - vertStart[z-1], 1);CHKERRQ(ierr);} 379f8716870SMatthew G. Knepley else {ierr = DMLabelSetValue(label, c - cellStart[z-1], 1);CHKERRQ(ierr);} 380f8716870SMatthew G. Knepley } 381c4cbe8f1SToby Isaac } else if (pointtype == CGNS_ENUMV(PointList)) { 382c4cbe8f1SToby Isaac CGNS_ENUMT(GridLocation_t) gridloc; 383f8716870SMatthew G. Knepley 384f8716870SMatthew G. Knepley /* List of points: Oh please, someone get the CGNS developers away from a computer. This is unconscionable. */ 385c261946bSLisandro Dalcin ierr = cg_goto(cgid, 1, "Zone_t", z, "BC_t", bc, "end");CHKERRCGNS(ierr); 386c261946bSLisandro Dalcin ierr = cg_gridlocation_read(&gridloc);CHKERRCGNS(ierr); 387f8716870SMatthew G. Knepley for (c = 0; c < npoints; ++c) { 388c4cbe8f1SToby Isaac if (gridloc == CGNS_ENUMV(Vertex)) {ierr = DMLabelSetValue(label, points[c] - vertStart[z-1], 1);CHKERRQ(ierr);} 389f8716870SMatthew G. Knepley else {ierr = DMLabelSetValue(label, points[c] - cellStart[z-1], 1);CHKERRQ(ierr);} 390f8716870SMatthew G. Knepley } 39198921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_SUP, "Unsupported point set type %d", (int) pointtype); 392f8716870SMatthew G. Knepley ierr = PetscFree2(points, normals);CHKERRQ(ierr); 393f8716870SMatthew G. Knepley } 394f8716870SMatthew G. Knepley } 395f8716870SMatthew G. Knepley ierr = PetscFree2(cellStart, vertStart);CHKERRQ(ierr); 396f8716870SMatthew G. Knepley } 3974e6cc61cSLisandro Dalcin ierr = DMGetNumLabels(*dm, &labelIdRange[1]);CHKERRQ(ierr); 3984e6cc61cSLisandro Dalcin ierr = MPI_Bcast(labelIdRange, 2, MPIU_INT, 0, comm);CHKERRMPI(ierr); 3994e6cc61cSLisandro Dalcin 4004e6cc61cSLisandro Dalcin /* Create BC labels at all processes */ 4014e6cc61cSLisandro Dalcin for (labelId = labelIdRange[0]; labelId < labelIdRange[1]; ++labelId) { 4024e6cc61cSLisandro Dalcin char *labelName = buffer; 4034e6cc61cSLisandro Dalcin size_t len = sizeof(buffer); 4044e6cc61cSLisandro Dalcin const char *locName; 4054e6cc61cSLisandro Dalcin 406dd400576SPatrick Sanan if (rank == 0) { 4074e6cc61cSLisandro Dalcin ierr = DMGetLabelByNum(*dm, labelId, &label);CHKERRQ(ierr); 4084e6cc61cSLisandro Dalcin ierr = PetscObjectGetName((PetscObject)label, &locName);CHKERRQ(ierr); 4094e6cc61cSLisandro Dalcin ierr = PetscStrncpy(labelName, locName, len);CHKERRQ(ierr); 4104e6cc61cSLisandro Dalcin } 4114e6cc61cSLisandro Dalcin ierr = MPI_Bcast(labelName, (PetscMPIInt)len, MPIU_INT, 0, comm);CHKERRMPI(ierr); 4124e6cc61cSLisandro Dalcin ierr = DMCreateLabel(*dm, labelName);CHKERRMPI(ierr); 4134e6cc61cSLisandro Dalcin } 414b458e8f1SJose E. Roman PetscFunctionReturn(0); 415e1b39ce3SMatthew G. Knepley #else 416e1b39ce3SMatthew G. Knepley SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir"); 417e1b39ce3SMatthew G. Knepley #endif 418e1b39ce3SMatthew G. Knepley } 419