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