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