1 #define PETSCDM_DLL 2 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 4 #if defined(PETSC_HAVE_CGNS) 5 #include <cgnslib.h> 6 #include <cgns_io.h> 7 #endif 8 9 #undef __FUNCT__ 10 #define __FUNCT__ "DMPlexCreateCGNS" 11 /*@ 12 DMPlexCreateCGNS - Create a DMPlex mesh from a CGNS file. 13 14 Collective on comm 15 16 Input Parameters: 17 + comm - The MPI communicator 18 . cgid - The CG id associated with a file and obtained using cg_open 19 - interpolate - Create faces and edges in the mesh 20 21 Output Parameter: 22 . dm - The DM object representing the mesh 23 24 Note: http://www.grc.nasa.gov/WWW/cgns/CGNS_docs_current/index.html 25 26 Level: beginner 27 28 .keywords: mesh,CGNS 29 .seealso: DMPlexCreate(), DMPlexCreateExodus() 30 @*/ 31 PetscErrorCode DMPlexCreateCGNS(MPI_Comm comm, PetscInt cgid, PetscBool interpolate, DM *dm) 32 { 33 #if defined(PETSC_HAVE_CGNS) 34 PetscMPIInt num_proc, rank; 35 PetscSection coordSection; 36 Vec coordinates; 37 PetscScalar *coords; 38 PetscInt coordSize, v; 39 PetscErrorCode ierr; 40 /* Read from file */ 41 char basename[CGIO_MAX_NAME_LENGTH+1]; 42 char buffer[CGIO_MAX_NAME_LENGTH+1]; 43 int dim = 0, physDim = 0, numVertices = 0, numCells = 0; 44 int nzones = 0; 45 #endif 46 47 PetscFunctionBegin; 48 #if defined(PETSC_HAVE_CGNS) 49 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 50 ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 51 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 52 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 53 /* Open CGNS II file and read basic informations on rank 0, then broadcast to all processors */ 54 if (!rank) { 55 int nbases, z; 56 57 ierr = cg_nbases(cgid, &nbases);CHKERRQ(ierr); 58 if (nbases > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single base, not %d\n",nbases); 59 ierr = cg_base_read(cgid, 1, basename, &dim, &physDim);CHKERRQ(ierr); 60 ierr = cg_nzones(cgid, 1, &nzones);CHKERRQ(ierr); 61 for (z = 1; z <= nzones; ++z) { 62 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 63 64 ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr); 65 numVertices += sizes[0]; 66 numCells += sizes[1]; 67 } 68 } 69 ierr = MPI_Bcast(basename, CGIO_MAX_NAME_LENGTH+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 70 ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 71 ierr = MPI_Bcast(&nzones, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 72 ierr = PetscObjectSetName((PetscObject) *dm, basename);CHKERRQ(ierr); 73 ierr = DMPlexSetDimension(*dm, dim);CHKERRQ(ierr); 74 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 75 76 /* Read zone information */ 77 if (!rank) { 78 int z, c, c_loc, v, v_loc; 79 80 /* Read the cell set connectivity table and build mesh topology 81 CGNS standard requires that cells in a zone be numbered sequentially and be pairwise disjoint. */ 82 /* First set sizes */ 83 for (z = 1, c = 0; z <= nzones; ++z) { 84 ZoneType_t zonetype; 85 int nsections; 86 ElementType_t cellType; 87 cgsize_t start, end; 88 int nbndry, parentFlag; 89 PetscInt numCorners; 90 91 ierr = cg_zone_type(cgid, 1, z, &zonetype);CHKERRQ(ierr); 92 if (zonetype == Structured) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_LIB,"Can only handle Unstructured zones for CGNS"); 93 ierr = cg_nsections(cgid, 1, z, &nsections);CHKERRQ(ierr); 94 if (nsections > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single section, not %d\n",nsections); 95 ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr); 96 /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 97 if (cellType == MIXED) { 98 cgsize_t elementDataSize, *elements; 99 PetscInt off; 100 101 ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr); 102 ierr = PetscMalloc(elementDataSize * sizeof(cgsize_t), &elements);CHKERRQ(ierr); 103 ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr); 104 for (c_loc = start, off = 0; c_loc < end; ++c_loc, ++c) { 105 cellType = elements[off]; 106 switch (cellType) { 107 case TRI_3: numCorners = 3;break; 108 case QUAD_4: numCorners = 4;break; 109 case TETRA_4: numCorners = 4;break; 110 case HEXA_8: numCorners = 8;break; 111 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 112 } 113 ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 114 off += numCorners+1; 115 } 116 ierr = PetscFree(elements);CHKERRQ(ierr); 117 } else { 118 switch (cellType) { 119 case TRI_3: numCorners = 3;break; 120 case QUAD_4: numCorners = 4;break; 121 case TETRA_4: numCorners = 4;break; 122 case HEXA_8: numCorners = 8;break; 123 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 124 } 125 for (c_loc = start; c_loc < end; ++c_loc, ++c) { 126 ierr = DMPlexSetConeSize(*dm, c, numCorners);CHKERRQ(ierr); 127 } 128 } 129 } 130 ierr = DMSetUp(*dm);CHKERRQ(ierr); 131 for (z = 1, c = 0; z <= nzones; ++z) { 132 ElementType_t cellType; 133 cgsize_t *elements, elementDataSize, start, end; 134 int nbndry, parentFlag; 135 PetscInt *cone, numc, numCorners, maxCorners = 27; 136 137 ierr = cg_section_read(cgid, 1, z, 1, buffer, &cellType, &start, &end, &nbndry, &parentFlag);CHKERRQ(ierr); 138 numc = end - start; 139 /* This alone is reason enough to bludgeon every single CGNDS developer, this must be what they describe as the "idiocy of crowds" */ 140 ierr = cg_ElementDataSize(cgid, 1, z, 1, &elementDataSize);CHKERRQ(ierr); 141 ierr = PetscMalloc2(elementDataSize,cgsize_t,&elements,maxCorners,PetscInt,&cone);CHKERRQ(ierr); 142 ierr = cg_elements_read(cgid, 1, z, 1, elements, NULL);CHKERRQ(ierr); 143 if (cellType == MIXED) { 144 /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */ 145 for (c_loc = 0, v = 0; c_loc < numc; ++c_loc, ++c) { 146 cellType = elements[v]; ++v; 147 switch (cellType) { 148 case TRI_3: numCorners = 3;break; 149 case QUAD_4: numCorners = 4;break; 150 case TETRA_4: numCorners = 4;break; 151 case HEXA_8: numCorners = 8;break; 152 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 153 } 154 for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 155 cone[v_loc] = elements[v]+numCells-1; 156 } 157 /* Tetrahedra are inverted */ 158 if (cellType == TETRA_4) { 159 PetscInt tmp = cone[0]; 160 cone[0] = cone[1]; 161 cone[1] = tmp; 162 } 163 /* Hexahedra are inverted */ 164 if (cellType == HEXA_8) { 165 PetscInt tmp = cone[1]; 166 cone[1] = cone[3]; 167 cone[3] = tmp; 168 } 169 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 170 ierr = DMPlexSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr); 171 } 172 } else { 173 switch (cellType) { 174 case TRI_3: numCorners = 3;break; 175 case QUAD_4: numCorners = 4;break; 176 case TETRA_4: numCorners = 4;break; 177 case HEXA_8: numCorners = 8;break; 178 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid cell type %d", (int) cellType); 179 } 180 181 /* CGNS uses Fortran-based indexing, sieve uses C-style and numbers cell first then vertices. */ 182 for (c_loc = 0, v = 0; c_loc < numc; ++c_loc, ++c) { 183 for (v_loc = 0; v_loc < numCorners; ++v_loc, ++v) { 184 cone[v_loc] = elements[v]+numCells-1; 185 } 186 /* Tetrahedra are inverted */ 187 if (cellType == TETRA_4) { 188 PetscInt tmp = cone[0]; 189 cone[0] = cone[1]; 190 cone[1] = tmp; 191 } 192 /* Hexahedra are inverted */ 193 if (cellType == HEXA_8) { 194 PetscInt tmp = cone[1]; 195 cone[1] = cone[3]; 196 cone[3] = tmp; 197 } 198 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 199 ierr = DMPlexSetLabelValue(*dm, "zone", c, z);CHKERRQ(ierr); 200 } 201 } 202 ierr = PetscFree2(elements,cone);CHKERRQ(ierr); 203 } 204 } 205 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 206 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 207 if (interpolate) { 208 DM idm; 209 210 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 211 /* Maintain zone label */ 212 { 213 DMLabel label; 214 215 ierr = DMPlexRemoveLabel(*dm, "zone", &label);CHKERRQ(ierr); 216 if (label) {ierr = DMPlexAddLabel(idm, label);CHKERRQ(ierr);} 217 } 218 ierr = DMDestroy(dm);CHKERRQ(ierr); 219 *dm = idm; 220 } 221 222 /* Read coordinates */ 223 ierr = DMPlexGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 224 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 225 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 226 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 227 for (v = numCells; v < numCells+numVertices; ++v) { 228 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 229 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 230 } 231 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 232 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 233 ierr = VecCreate(comm, &coordinates);CHKERRQ(ierr); 234 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 235 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 236 ierr = VecSetFromOptions(coordinates);CHKERRQ(ierr); 237 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 238 if (!rank) { 239 PetscInt off = 0; 240 float *x[3]; 241 int z, c, d; 242 243 ierr = PetscMalloc3(numVertices,float,&x[0],numVertices,float,&x[1],numVertices,float,&x[2]);CHKERRQ(ierr); 244 for (z = 1, c = 0; z <= nzones; ++z) { 245 DataType_t datatype; 246 cgsize_t sizes[3]; /* Number of vertices, number of cells, number of boundary vertices */ 247 cgsize_t range_min[3] = {1, 1, 1}; 248 cgsize_t range_max[3] = {1, 1, 1}; 249 int ngrids, ncoords; 250 251 252 ierr = cg_zone_read(cgid, 1, z, buffer, sizes);CHKERRQ(ierr); 253 range_max[0] = sizes[0]; 254 ierr = cg_ngrids(cgid, 1, z, &ngrids);CHKERRQ(ierr); 255 if (ngrids > 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a single grid, not %d\n",ngrids); 256 ierr = cg_ncoords(cgid, 1, z, &ncoords);CHKERRQ(ierr); 257 if (ncoords != dim) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"CGNS file must have a coordinate array for each dimension, not %d\n",ncoords); 258 for (d = 0; d < dim; ++d) { 259 ierr = cg_coord_info(cgid, 1, z, 1, &datatype, buffer);CHKERRQ(ierr); 260 ierr = cg_coord_read(cgid, 1, z, buffer, RealSingle, range_min, range_max, x[d]);CHKERRQ(ierr); 261 } 262 if (dim > 0) { 263 for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+0] = x[0][v]; 264 } 265 if (dim > 1) { 266 for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+1] = x[1][v]; 267 } 268 if (dim > 2) { 269 for (v = 0; v < sizes[0]; ++v) coords[(v+off)*dim+2] = x[2][v]; 270 } 271 off += sizes[0]; 272 } 273 ierr = PetscFree3(x[0],x[1],x[2]);CHKERRQ(ierr); 274 } 275 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 276 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 277 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 278 #else 279 SETERRQ(comm, PETSC_ERR_SUP, "This method requires CGNS support. Reconfigure using --with-cgns-dir"); 280 #endif 281 PetscFunctionReturn(0); 282 } 283