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