1 #define PETSCDM_DLL 2 #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3 4 #if defined(PETSC_HAVE_EXODUSII) 5 #include <netcdf.h> 6 #include <exodusII.h> 7 #endif 8 9 #if defined(PETSC_HAVE_EXODUSII) 10 /* 11 EXOGetVarIndex - Locate a result in an exodus file based on its name 12 13 Not collective 14 15 Input Parameters: 16 + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 17 . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 18 - name - the name of the result 19 20 Output Parameters: 21 . varIndex - the location in the exodus file of the result 22 23 Notes: 24 The exodus variable index is obtained by comparing name and the 25 names of zonal variables declared in the exodus file. For instance if name is "V" 26 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 27 amongst all variables of type obj_type. 28 29 Level: beginner 30 31 .keywords: mesh,ExodusII 32 .seealso: EXOGetVarIndex(),DMView_Exodus(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO() 33 */ 34 static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 35 { 36 int num_vars, i, j; 37 char ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1]; 38 const int num_suffix = 5; 39 char *suffix[5]; 40 PetscBool flg; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 suffix[0] = (char *) ""; 45 suffix[1] = (char *) "_X"; 46 suffix[2] = (char *) "_XX"; 47 suffix[3] = (char *) "_1"; 48 suffix[4] = (char *) "_11"; 49 50 *varIndex = 0; 51 PetscStackCallStandard(ex_get_variable_param,(exoid, obj_type, &num_vars)); 52 for (i = 0; i < num_vars; ++i) { 53 PetscStackCallStandard(ex_get_variable_name,(exoid, obj_type, i+1, var_name)); 54 for (j = 0; j < num_suffix; ++j){ 55 ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr); 56 ierr = PetscStrncat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr); 57 ierr = PetscStrcasecmp(ext_name, var_name, &flg); 58 if (flg) { 59 *varIndex = i+1; 60 PetscFunctionReturn(0); 61 } 62 } 63 } 64 SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name); 65 PetscFunctionReturn(-1); 66 } 67 68 /* 69 DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format 70 71 Collective on dm 72 73 Input Parameters: 74 + dm - The dm to be written 75 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 76 - degree - the degree of the interpolation space 77 78 Notes: 79 Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 80 consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 81 82 If the dm has been distributed and exoid points to different files on each MPI rank, only the "local" part of the DM 83 (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will 84 probably be corrupted. 85 86 DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 87 on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 88 It should be extended to use PetscFE objects. 89 90 This function will only handle TRI, TET, QUAD and HEX cells. 91 Level: beginner 92 93 .keywords: mesh, ExodusII 94 .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 95 */ 96 PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree) 97 { 98 enum ElemType {TRI, QUAD, TET, HEX}; 99 MPI_Comm comm; 100 /* Connectivity Variables */ 101 PetscInt cellsNotInConnectivity; 102 /* Cell Sets */ 103 DMLabel csLabel; 104 IS csIS; 105 const PetscInt *csIdx; 106 PetscInt num_cs, cs; 107 enum ElemType *type; 108 /* Coordinate Variables */ 109 DM cdm; 110 PetscSection section; 111 Vec coord; 112 PetscInt **nodes; 113 PetscInt depth, d, dim, skipCells = 0; 114 PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 115 PetscInt num_vs, num_fs; 116 PetscMPIInt rank, size; 117 const char *dmName; 118 PetscErrorCode ierr; 119 120 PetscFunctionBegin; 121 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 122 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 123 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 124 /* --- Get DM info --- */ 125 ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr); 126 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 127 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 128 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 129 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 130 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 131 ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr); 132 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 133 numCells = cEnd - cStart; 134 numEdges = eEnd - eStart; 135 numVertices = vEnd - vStart; 136 if (depth == 3) {numFaces = fEnd - fStart;} 137 else {numFaces = 0;} 138 ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr); 139 ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr); 140 ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr); 141 ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 142 ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 143 if (num_cs > 0) { 144 ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr); 145 ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr); 146 ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr); 147 } 148 ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr); 149 /* Set element type for each block and compute total number of nodes */ 150 ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr); 151 numNodes = numVertices; 152 if (degree == 2) {numNodes += numEdges;} 153 cellsNotInConnectivity = numCells; 154 for (cs = 0; cs < num_cs; ++cs) { 155 IS stratumIS; 156 const PetscInt *cells; 157 PetscScalar *xyz = NULL; 158 PetscInt csSize, closureSize; 159 PetscInt nodesTriP1[4] = {3,0,0,0}; 160 PetscInt nodesTriP2[4] = {3,3,0,0}; 161 PetscInt nodesQuadP1[4] = {4,0,0,0}; 162 PetscInt nodesQuadP2[4] = {4,4,0,1}; 163 PetscInt nodesTetP1[4] = {4,0,0,0}; 164 PetscInt nodesTetP2[4] = {4,6,0,0}; 165 PetscInt nodesHexP1[4] = {8,0,0,0}; 166 PetscInt nodesHexP2[4] = {8,12,6,1}; 167 168 ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 169 ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 170 ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 171 ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 172 switch (dim) { 173 case 2: 174 if (closureSize == 3*dim) {type[cs] = TRI;} 175 else if (closureSize == 4*dim) {type[cs] = QUAD;} 176 else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 177 break; 178 case 3: 179 if (closureSize == 4*dim) {type[cs] = TET;} 180 else if (closureSize == 8*dim) {type[cs] = HEX;} 181 else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 182 break; 183 default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim); 184 } 185 if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;} 186 if ((degree == 2) && (type[cs] == HEX)) {numNodes += csSize; numNodes += numFaces;} 187 ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 188 /* Set nodes and Element type */ 189 if (type[cs] == TRI) { 190 if (degree == 1) nodes[cs] = nodesTriP1; 191 else if (degree == 2) nodes[cs] = nodesTriP2; 192 } else if (type[cs] == QUAD) { 193 if (degree == 1) nodes[cs] = nodesQuadP1; 194 else if (degree == 2) nodes[cs] = nodesQuadP2; 195 } else if (type[cs] == TET) { 196 if (degree == 1) nodes[cs] = nodesTetP1; 197 else if (degree == 2) nodes[cs] = nodesTetP2; 198 } else if (type[cs] == HEX) { 199 if (degree == 1) nodes[cs] = nodesHexP1; 200 else if (degree == 2) nodes[cs] = nodesHexP2; 201 } 202 /* Compute the number of cells not in the connectivity table */ 203 cellsNotInConnectivity -= nodes[cs][3]*csSize; 204 205 ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 206 ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 207 } 208 if (num_cs > 0) {PetscStackCallStandard(ex_put_init,(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs));} 209 /* --- Connectivity --- */ 210 for (cs = 0; cs < num_cs; ++cs) { 211 IS stratumIS; 212 const PetscInt *cells; 213 PetscInt *connect, off = 0; 214 PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 215 PetscInt csSize, c, connectSize, closureSize; 216 char *elem_type = NULL; 217 char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 218 char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 219 char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 220 char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 221 222 ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 223 ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 224 ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 225 /* Set Element type */ 226 if (type[cs] == TRI) { 227 if (degree == 1) elem_type = elem_type_tri3; 228 else if (degree == 2) elem_type = elem_type_tri6; 229 } else if (type[cs] == QUAD) { 230 if (degree == 1) elem_type = elem_type_quad4; 231 else if (degree == 2) elem_type = elem_type_quad9; 232 } else if (type[cs] == TET) { 233 if (degree == 1) elem_type = elem_type_tet4; 234 else if (degree == 2) elem_type = elem_type_tet10; 235 } else if (type[cs] == HEX) { 236 if (degree == 1) elem_type = elem_type_hex8; 237 else if (degree == 2) elem_type = elem_type_hex27; 238 } 239 connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 240 ierr = PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect);CHKERRQ(ierr); 241 PetscStackCallStandard(ex_put_block,(exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1)); 242 /* Find number of vertices, edges, and faces in the closure */ 243 verticesInClosure = nodes[cs][0]; 244 if (depth > 1) { 245 if (dim == 2) { 246 ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr); 247 } else if (dim == 3) { 248 PetscInt *closure = NULL; 249 250 ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr); 251 ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 252 edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 253 ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 254 } 255 } 256 /* Get connectivity for each cell */ 257 for (c = 0; c < csSize; ++c) { 258 PetscInt *closure = NULL; 259 PetscInt temp, i; 260 261 ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 262 for (i = 0; i < connectSize; ++i) { 263 if (i < nodes[cs][0]) {/* Vertices */ 264 connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 265 connect[i+off] -= cellsNotInConnectivity; 266 } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 267 connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 268 if (nodes[cs][2] == 0) connect[i+off] -= numFaces; 269 connect[i+off] -= cellsNotInConnectivity; 270 } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */ 271 connect[i+off] = closure[0] + 1; 272 connect[i+off] -= skipCells; 273 } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */ 274 connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 275 connect[i+off] -= cellsNotInConnectivity; 276 } else { 277 connect[i+off] = -1; 278 } 279 } 280 /* Tetrahedra are inverted */ 281 if (type[cs] == TET) { 282 temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp; 283 if (degree == 2) { 284 temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp; 285 temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp; 286 } 287 } 288 /* Hexahedra are inverted */ 289 if (type[cs] == HEX) { 290 temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp; 291 if (degree == 2) { 292 temp = connect[8+off]; connect[8+off] = connect[11+off]; connect[11+off] = temp; 293 temp = connect[9+off]; connect[9+off] = connect[10+off]; connect[10+off] = temp; 294 temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp; 295 temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp; 296 297 temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp; 298 temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp; 299 temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp; 300 temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp; 301 302 temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp; 303 temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp; 304 temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp; 305 } 306 } 307 off += connectSize; 308 ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 309 } 310 PetscStackCallStandard(ex_put_conn,(exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0)); 311 skipCells += (nodes[cs][3] == 0)*csSize; 312 ierr = PetscFree(connect);CHKERRQ(ierr); 313 ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 314 ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 315 } 316 ierr = PetscFree(type);CHKERRQ(ierr); 317 /* --- Coordinates --- */ 318 ierr = PetscSectionCreate(comm, §ion);CHKERRQ(ierr); 319 ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 320 for (d = 0; d < depth; ++d) { 321 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 322 for (p = pStart; p < pEnd; ++p) { 323 ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr); 324 } 325 } 326 for (cs = 0; cs < num_cs; ++cs) { 327 IS stratumIS; 328 const PetscInt *cells; 329 PetscInt csSize, c; 330 331 ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 332 ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 333 ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 334 for (c = 0; c < csSize; ++c) { 335 ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr); 336 } 337 ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 338 ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 339 } 340 if (num_cs > 0) { 341 ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr); 342 ierr = ISDestroy(&csIS);CHKERRQ(ierr); 343 } 344 ierr = PetscFree(nodes);CHKERRQ(ierr); 345 ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 346 if (numNodes > 0) { 347 const char *coordNames[3] = {"x", "y", "z"}; 348 PetscScalar *coords, *closure; 349 PetscReal *cval; 350 PetscInt hasDof, n = 0; 351 352 /* There can't be more than 24 values in the closure of a point for the coord section */ 353 ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr); 354 ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 355 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 356 for (p = pStart; p < pEnd; ++p) { 357 ierr = PetscSectionGetDof(section, p, &hasDof); 358 if (hasDof) { 359 PetscInt closureSize = 24, j; 360 361 ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr); 362 for (d = 0; d < dim; ++d) { 363 cval[d] = 0.0; 364 for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 365 coords[d*numNodes+n] = cval[d] * dim / closureSize; 366 } 367 ++n; 368 } 369 } 370 PetscStackCallStandard(ex_put_coord,(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes])); 371 ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr); 372 PetscStackCallStandard(ex_put_coord_names,(exoid, (char **) coordNames)); 373 } 374 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 /* 379 VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 380 381 Collective on v 382 383 Input Parameters: 384 + v - The vector to be written 385 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 386 - step - the time step to write at (exodus steps are numbered starting from 1) 387 388 Notes: 389 The exodus result nodal variable index is obtained by comparing the Vec name and the 390 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 391 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 392 amongst all nodal variables. 393 394 Level: beginner 395 396 .keywords: mesh,ExodusII 397 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO() 398 @*/ 399 PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 400 { 401 MPI_Comm comm; 402 PetscMPIInt size; 403 DM dm; 404 Vec vNatural, vComp; 405 const PetscReal *varray; 406 const char *vecname; 407 PetscInt xs, xe, bs; 408 PetscBool useNatural; 409 int offset; 410 PetscErrorCode ierr; 411 412 PetscFunctionBegin; 413 ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 414 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 415 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 416 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 417 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 418 if (useNatural) { 419 ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 420 ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 421 ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 422 } else { 423 vNatural = v; 424 } 425 /* Get the location of the variable in the exodus file */ 426 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 427 ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 428 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname); 429 /* Write local chunk of the result in the exodus file 430 exodus stores each component of a vector-valued field as a separate variable. 431 We assume that they are stored sequentially */ 432 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 433 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 434 if (bs == 1) { 435 ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 436 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 437 ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 438 } else { 439 IS compIS; 440 PetscInt c; 441 442 ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 443 for (c = 0; c < bs; ++c) { 444 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 445 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 446 ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 447 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 448 ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 449 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 450 } 451 ierr = ISDestroy(&compIS);CHKERRQ(ierr); 452 } 453 if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);} 454 PetscFunctionReturn(0); 455 } 456 457 /* 458 VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 459 460 Collective on v 461 462 Input Parameters: 463 + v - The vector to be written 464 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 465 - step - the time step to read at (exodus steps are numbered starting from 1) 466 467 Notes: 468 The exodus result nodal variable index is obtained by comparing the Vec name and the 469 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 470 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 471 amongst all nodal variables. 472 473 Level: beginner 474 475 .keywords: mesh, ExodusII 476 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 477 */ 478 PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 479 { 480 MPI_Comm comm; 481 PetscMPIInt size; 482 DM dm; 483 Vec vNatural, vComp; 484 PetscReal *varray; 485 const char *vecname; 486 PetscInt xs, xe, bs; 487 PetscBool useNatural; 488 int offset; 489 PetscErrorCode ierr; 490 491 PetscFunctionBegin; 492 ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 493 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 494 ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 495 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 496 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 497 if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 498 else {vNatural = v;} 499 /* Get the location of the variable in the exodus file */ 500 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 501 ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 502 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname); 503 /* Read local chunk from the file */ 504 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 505 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 506 if (bs == 1) { 507 ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 508 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 509 ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 510 } else { 511 IS compIS; 512 PetscInt c; 513 514 ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 515 for (c = 0; c < bs; ++c) { 516 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 517 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 518 ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 519 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 520 ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 521 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 522 } 523 ierr = ISDestroy(&compIS);CHKERRQ(ierr); 524 } 525 if (useNatural) { 526 ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 527 ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 528 ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 529 } 530 PetscFunctionReturn(0); 531 } 532 533 /* 534 VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 535 536 Collective on v 537 538 Input Parameters: 539 + v - The vector to be written 540 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 541 - step - the time step to write at (exodus steps are numbered starting from 1) 542 543 Notes: 544 The exodus result zonal variable index is obtained by comparing the Vec name and the 545 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 546 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 547 amongst all zonal variables. 548 549 Level: beginner 550 551 .keywords: mesh,ExodusII 552 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal() 553 */ 554 PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 555 { 556 MPI_Comm comm; 557 PetscMPIInt size; 558 DM dm; 559 Vec vNatural, vComp; 560 const PetscReal *varray; 561 const char *vecname; 562 PetscInt xs, xe, bs; 563 PetscBool useNatural; 564 int offset; 565 IS compIS; 566 PetscInt *csSize, *csID; 567 PetscInt numCS, set, csxs = 0; 568 PetscErrorCode ierr; 569 570 PetscFunctionBegin; 571 ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr); 572 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 573 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 574 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 575 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 576 if (useNatural) { 577 ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 578 ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 579 ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 580 } else { 581 vNatural = v; 582 } 583 /* Get the location of the variable in the exodus file */ 584 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 585 ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 586 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 587 /* Write local chunk of the result in the exodus file 588 exodus stores each component of a vector-valued field as a separate variable. 589 We assume that they are stored sequentially 590 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 591 but once the vector has been reordered to natural size, we cannot use the label informations 592 to figure out what to save where. */ 593 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 594 ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 595 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 596 for (set = 0; set < numCS; ++set) { 597 ex_block block; 598 599 block.id = csID[set]; 600 block.type = EX_ELEM_BLOCK; 601 PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 602 csSize[set] = block.num_entry; 603 } 604 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 605 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 606 if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 607 for (set = 0; set < numCS; set++) { 608 PetscInt csLocalSize, c; 609 610 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 611 local slice of zonal values: xs/bs,xm/bs-1 612 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 613 csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 614 if (bs == 1) { 615 ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 616 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 617 ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 618 } else { 619 for (c = 0; c < bs; ++c) { 620 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 621 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 622 ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 623 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)])); 624 ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 625 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 626 } 627 } 628 csxs += csSize[set]; 629 } 630 ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 631 if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 632 if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 633 PetscFunctionReturn(0); 634 } 635 636 /* 637 VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 638 639 Collective on v 640 641 Input Parameters: 642 + v - The vector to be written 643 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 644 - step - the time step to read at (exodus steps are numbered starting from 1) 645 646 Notes: 647 The exodus result zonal variable index is obtained by comparing the Vec name and the 648 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 649 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 650 amongst all zonal variables. 651 652 Level: beginner 653 654 .keywords: mesh,ExodusII 655 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal() 656 */ 657 PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 658 { 659 MPI_Comm comm; 660 PetscMPIInt size; 661 DM dm; 662 Vec vNatural, vComp; 663 PetscReal *varray; 664 const char *vecname; 665 PetscInt xs, xe, bs; 666 PetscBool useNatural; 667 int offset; 668 IS compIS; 669 PetscInt *csSize, *csID; 670 PetscInt numCS, set, csxs = 0; 671 PetscErrorCode ierr; 672 673 PetscFunctionBegin; 674 ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 675 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 676 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 677 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 678 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 679 if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 680 else {vNatural = v;} 681 /* Get the location of the variable in the exodus file */ 682 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 683 ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 684 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 685 /* Read local chunk of the result in the exodus file 686 exodus stores each component of a vector-valued field as a separate variable. 687 We assume that they are stored sequentially 688 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 689 but once the vector has been reordered to natural size, we cannot use the label informations 690 to figure out what to save where. */ 691 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 692 ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 693 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 694 for (set = 0; set < numCS; ++set) { 695 ex_block block; 696 697 block.id = csID[set]; 698 block.type = EX_ELEM_BLOCK; 699 PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 700 csSize[set] = block.num_entry; 701 } 702 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 703 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 704 if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 705 for (set = 0; set < numCS; ++set) { 706 PetscInt csLocalSize, c; 707 708 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 709 local slice of zonal values: xs/bs,xm/bs-1 710 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 711 csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 712 if (bs == 1) { 713 ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 714 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 715 ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 716 } else { 717 for (c = 0; c < bs; ++c) { 718 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 719 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 720 ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 721 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)])); 722 ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 723 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 724 } 725 } 726 csxs += csSize[set]; 727 } 728 ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 729 if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 730 if (useNatural) { 731 ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 732 ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 733 ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 734 } 735 PetscFunctionReturn(0); 736 } 737 #endif 738 739 /*@C 740 DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 741 742 Collective on comm 743 744 Input Parameters: 745 + comm - The MPI communicator 746 . filename - The name of the ExodusII file 747 - interpolate - Create faces and edges in the mesh 748 749 Output Parameter: 750 . dm - The DM object representing the mesh 751 752 Level: beginner 753 754 .keywords: mesh,ExodusII 755 .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 756 @*/ 757 PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 758 { 759 PetscMPIInt rank; 760 PetscErrorCode ierr; 761 #if defined(PETSC_HAVE_EXODUSII) 762 int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 763 float version; 764 #endif 765 766 PetscFunctionBegin; 767 PetscValidCharPointer(filename, 2); 768 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 769 #if defined(PETSC_HAVE_EXODUSII) 770 if (!rank) { 771 exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 772 if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 773 } 774 ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 775 if (!rank) {PetscStackCallStandard(ex_close,(exoid));} 776 #else 777 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 778 #endif 779 PetscFunctionReturn(0); 780 } 781 782 /*@ 783 DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 784 785 Collective on comm 786 787 Input Parameters: 788 + comm - The MPI communicator 789 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 790 - interpolate - Create faces and edges in the mesh 791 792 Output Parameter: 793 . dm - The DM object representing the mesh 794 795 Level: beginner 796 797 .keywords: mesh,ExodusII 798 .seealso: DMPLEX, DMCreate() 799 @*/ 800 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 801 { 802 #if defined(PETSC_HAVE_EXODUSII) 803 PetscMPIInt num_proc, rank; 804 PetscSection coordSection; 805 Vec coordinates; 806 PetscScalar *coords; 807 PetscInt coordSize, v; 808 PetscErrorCode ierr; 809 /* Read from ex_get_init() */ 810 char title[PETSC_MAX_PATH_LEN+1]; 811 int dim = 0, numVertices = 0, numCells = 0; 812 int num_cs = 0, num_vs = 0, num_fs = 0; 813 #endif 814 815 PetscFunctionBegin; 816 #if defined(PETSC_HAVE_EXODUSII) 817 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 818 ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 819 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 820 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 821 /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 822 if (!rank) { 823 ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 824 PetscStackCallStandard(ex_get_init,(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs)); 825 if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 826 } 827 ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 828 ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 829 ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 830 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 831 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 832 833 /* Read cell sets information */ 834 if (!rank) { 835 PetscInt *cone; 836 int c, cs, c_loc, v, v_loc; 837 /* Read from ex_get_elem_blk_ids() */ 838 int *cs_id; 839 /* Read from ex_get_elem_block() */ 840 char buffer[PETSC_MAX_PATH_LEN+1]; 841 int num_cell_in_set, num_vertex_per_cell, num_attr; 842 /* Read from ex_get_elem_conn() */ 843 int *cs_connect; 844 845 /* Get cell sets IDs */ 846 ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr); 847 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id)); 848 /* Read the cell set connectivity table and build mesh topology 849 EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 850 /* First set sizes */ 851 for (cs = 0, c = 0; cs < num_cs; cs++) { 852 PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr)); 853 for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 854 ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 855 } 856 } 857 ierr = DMSetUp(*dm);CHKERRQ(ierr); 858 for (cs = 0, c = 0; cs < num_cs; cs++) { 859 PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr)); 860 ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 861 PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL)); 862 /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 863 for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 864 for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 865 cone[v_loc] = cs_connect[v]+numCells-1; 866 } 867 if (dim == 3) { 868 /* Tetrahedra are inverted */ 869 if (num_vertex_per_cell == 4) { 870 PetscInt tmp = cone[0]; 871 cone[0] = cone[1]; 872 cone[1] = tmp; 873 } 874 /* Hexahedra are inverted */ 875 if (num_vertex_per_cell == 8) { 876 PetscInt tmp = cone[1]; 877 cone[1] = cone[3]; 878 cone[3] = tmp; 879 } 880 } 881 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 882 ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 883 } 884 ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 885 } 886 ierr = PetscFree(cs_id);CHKERRQ(ierr); 887 } 888 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 889 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 890 if (interpolate) { 891 DM idm; 892 893 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 894 ierr = DMDestroy(dm);CHKERRQ(ierr); 895 *dm = idm; 896 } 897 898 /* Create vertex set label */ 899 if (!rank && (num_vs > 0)) { 900 int vs, v; 901 /* Read from ex_get_node_set_ids() */ 902 int *vs_id; 903 /* Read from ex_get_node_set_param() */ 904 int num_vertex_in_set; 905 /* Read from ex_get_node_set() */ 906 int *vs_vertex_list; 907 908 /* Get vertex set ids */ 909 ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 910 PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id)); 911 for (vs = 0; vs < num_vs; ++vs) { 912 PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL)); 913 ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 914 PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL)); 915 for (v = 0; v < num_vertex_in_set; ++v) { 916 ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 917 } 918 ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 919 } 920 ierr = PetscFree(vs_id);CHKERRQ(ierr); 921 } 922 /* Read coordinates */ 923 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 924 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 925 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 926 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 927 for (v = numCells; v < numCells+numVertices; ++v) { 928 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 929 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 930 } 931 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 932 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 933 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 934 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 935 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 936 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 937 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 938 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 939 if (!rank) { 940 PetscReal *x, *y, *z; 941 942 ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 943 PetscStackCallStandard(ex_get_coord,(exoid, x, y, z)); 944 if (dim > 0) { 945 for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 946 } 947 if (dim > 1) { 948 for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 949 } 950 if (dim > 2) { 951 for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 952 } 953 ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 954 } 955 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 956 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 957 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 958 959 /* Create side set label */ 960 if (!rank && interpolate && (num_fs > 0)) { 961 int fs, f, voff; 962 /* Read from ex_get_side_set_ids() */ 963 int *fs_id; 964 /* Read from ex_get_side_set_param() */ 965 int num_side_in_set; 966 /* Read from ex_get_side_set_node_list() */ 967 int *fs_vertex_count_list, *fs_vertex_list; 968 /* Read side set labels */ 969 char fs_name[MAX_STR_LENGTH+1]; 970 971 /* Get side set ids */ 972 ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 973 PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id)); 974 for (fs = 0; fs < num_fs; ++fs) { 975 PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL)); 976 ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 977 PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list)); 978 /* Get the specific name associated with this side set ID. */ 979 int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 980 for (f = 0, voff = 0; f < num_side_in_set; ++f) { 981 const PetscInt *faces = NULL; 982 PetscInt faceSize = fs_vertex_count_list[f], numFaces; 983 PetscInt faceVertices[4], v; 984 985 if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 986 for (v = 0; v < faceSize; ++v, ++voff) { 987 faceVertices[v] = fs_vertex_list[voff]+numCells-1; 988 } 989 ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 990 if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 991 ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 992 /* Only add the label if one has been detected for this side set. */ 993 if (!fs_name_err) { 994 ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 995 } 996 ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 997 } 998 ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 999 } 1000 ierr = PetscFree(fs_id);CHKERRQ(ierr); 1001 } 1002 #else 1003 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1004 #endif 1005 PetscFunctionReturn(0); 1006 } 1007