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; 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), &connect);CHKERRQ(ierr); 241 PetscStackCallStandard(ex_put_block,(exoid, EX_ELEM_BLOCK, 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] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 265 connect[i] -= cellsNotInConnectivity; 266 } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 267 connect[i] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 268 if (nodes[cs][2] == 0) connect[i] -= numFaces; 269 connect[i] -= cellsNotInConnectivity; 270 } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */ 271 connect[i] = closure[0] + 1; 272 connect[i] -= skipCells; 273 } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */ 274 connect[i] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 275 connect[i] -= cellsNotInConnectivity; 276 } else { 277 connect[i] = -1; 278 } 279 } 280 /* Tetrahedra are inverted */ 281 if (type[cs] == TET) { 282 temp = connect[0]; connect[0] = connect[1]; connect[1] = temp; 283 if (degree == 2) { 284 temp = connect[5]; connect[5] = connect[6]; connect[6] = temp; 285 temp = connect[7]; connect[7] = connect[8]; connect[8] = temp; 286 } 287 } 288 /* Hexahedra are inverted */ 289 if (type[cs] == HEX) { 290 temp = connect[1]; connect[1] = connect[3]; connect[3] = temp; 291 if (degree == 2) { 292 temp = connect[8]; connect[8] = connect[11]; connect[11] = temp; 293 temp = connect[9]; connect[9] = connect[10]; connect[10] = temp; 294 temp = connect[16]; connect[16] = connect[17]; connect[17] = temp; 295 temp = connect[18]; connect[18] = connect[19]; connect[19] = temp; 296 297 temp = connect[12]; connect[12] = connect[16]; connect[16] = temp; 298 temp = connect[13]; connect[13] = connect[17]; connect[17] = temp; 299 temp = connect[14]; connect[14] = connect[18]; connect[18] = temp; 300 temp = connect[15]; connect[15] = connect[19]; connect[19] = temp; 301 302 temp = connect[23]; connect[23] = connect[26]; connect[26] = temp; 303 temp = connect[24]; connect[24] = connect[25]; connect[25] = temp; 304 temp = connect[25]; connect[25] = connect[26]; connect[26] = temp; 305 } 306 } 307 PetscStackCallStandard(ex_put_partial_conn,(exoid, EX_ELEM_BLOCK, cs, c+1, 1, connect, NULL, NULL)); 308 ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 309 } 310 skipCells += (nodes[cs][3] == 0)*csSize; 311 ierr = PetscFree(connect);CHKERRQ(ierr); 312 ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 313 ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 314 } 315 ierr = PetscFree(type);CHKERRQ(ierr); 316 /* --- Coordinates --- */ 317 ierr = PetscSectionCreate(comm, §ion);CHKERRQ(ierr); 318 ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 319 for (d = 0; d < depth; ++d) { 320 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 321 for (p = pStart; p < pEnd; ++p) { 322 ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr); 323 } 324 } 325 for (cs = 0; cs < num_cs; ++cs) { 326 IS stratumIS; 327 const PetscInt *cells; 328 PetscInt csSize, c; 329 330 ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 331 ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 332 ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 333 for (c = 0; c < csSize; ++c) { 334 ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr); 335 } 336 ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 337 ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 338 } 339 if (num_cs > 0) { 340 ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr); 341 ierr = ISDestroy(&csIS);CHKERRQ(ierr); 342 } 343 ierr = PetscFree(nodes);CHKERRQ(ierr); 344 ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 345 if (numNodes > 0) { 346 const char *coordNames[3] = {"x", "y", "z"}; 347 PetscScalar *coords, *closure; 348 PetscReal *cval; 349 PetscInt hasDof, n = 0; 350 351 /* There can't be more than 24 values in the closure of a point for the coord section */ 352 ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr); 353 ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 354 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 355 for (p = pStart; p < pEnd; ++p) { 356 ierr = PetscSectionGetDof(section, p, &hasDof); 357 if (hasDof) { 358 PetscInt closureSize = 24, j; 359 360 ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr); 361 for (d = 0; d < dim; ++d) { 362 cval[d] = 0.0; 363 for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 364 coords[d*numNodes+n] = cval[d] * dim / closureSize; 365 } 366 ++n; 367 } 368 } 369 PetscStackCallStandard(ex_put_coord,(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes])); 370 ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr); 371 PetscStackCallStandard(ex_put_coord_names,(exoid, (char **) coordNames)); 372 } 373 ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 374 PetscFunctionReturn(0); 375 } 376 377 /* 378 VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 379 380 Collective on v 381 382 Input Parameters: 383 + v - The vector to be written 384 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 385 - step - the time step to write at (exodus steps are numbered starting from 1) 386 387 Notes: 388 The exodus result nodal variable index is obtained by comparing the Vec name and the 389 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 390 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 391 amongst all nodal variables. 392 393 Level: beginner 394 395 .keywords: mesh,ExodusII 396 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO() 397 @*/ 398 PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 399 { 400 MPI_Comm comm; 401 PetscMPIInt size; 402 DM dm; 403 Vec vNatural, vComp; 404 const PetscReal *varray; 405 const char *vecname; 406 PetscInt xs, xe, bs; 407 PetscBool useNatural; 408 int offset; 409 PetscErrorCode ierr; 410 411 PetscFunctionBegin; 412 ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 413 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 414 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 415 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 416 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 417 if (useNatural) { 418 ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 419 ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 420 ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 421 } else { 422 vNatural = v; 423 } 424 /* Get the location of the variable in the exodus file */ 425 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 426 ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 427 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname); 428 /* Write local chunk of the result in the exodus file 429 exodus stores each component of a vector-valued field as a separate variable. 430 We assume that they are stored sequentially */ 431 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 432 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 433 if (bs == 1) { 434 ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 435 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 436 ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 437 } else { 438 IS compIS; 439 PetscInt c; 440 441 ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 442 for (c = 0; c < bs; ++c) { 443 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 444 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 445 ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 446 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 447 ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 448 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 449 } 450 ierr = ISDestroy(&compIS);CHKERRQ(ierr); 451 } 452 if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);} 453 PetscFunctionReturn(0); 454 } 455 456 /* 457 VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 458 459 Collective on v 460 461 Input Parameters: 462 + v - The vector to be written 463 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 464 - step - the time step to read at (exodus steps are numbered starting from 1) 465 466 Notes: 467 The exodus result nodal variable index is obtained by comparing the Vec name and the 468 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 469 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 470 amongst all nodal variables. 471 472 Level: beginner 473 474 .keywords: mesh, ExodusII 475 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 476 */ 477 PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 478 { 479 MPI_Comm comm; 480 PetscMPIInt size; 481 DM dm; 482 Vec vNatural, vComp; 483 PetscReal *varray; 484 const char *vecname; 485 PetscInt xs, xe, bs; 486 PetscBool useNatural; 487 int offset; 488 PetscErrorCode ierr; 489 490 PetscFunctionBegin; 491 ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 492 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 493 ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 494 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 495 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 496 if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 497 else {vNatural = v;} 498 /* Get the location of the variable in the exodus file */ 499 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 500 ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 501 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname); 502 /* Read local chunk from the file */ 503 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 504 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 505 if (bs == 1) { 506 ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 507 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 508 ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 509 } else { 510 IS compIS; 511 PetscInt c; 512 513 ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 514 for (c = 0; c < bs; ++c) { 515 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 516 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 517 ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 518 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 519 ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 520 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 521 } 522 ierr = ISDestroy(&compIS);CHKERRQ(ierr); 523 } 524 if (useNatural) { 525 ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 526 ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 527 ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 528 } 529 PetscFunctionReturn(0); 530 } 531 532 /* 533 VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 534 535 Collective on v 536 537 Input Parameters: 538 + v - The vector to be written 539 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 540 - step - the time step to write at (exodus steps are numbered starting from 1) 541 542 Notes: 543 The exodus result zonal variable index is obtained by comparing the Vec name and the 544 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 545 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 546 amongst all zonal variables. 547 548 Level: beginner 549 550 .keywords: mesh,ExodusII 551 .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal() 552 */ 553 PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 554 { 555 MPI_Comm comm; 556 PetscMPIInt size; 557 DM dm; 558 Vec vNatural, vComp; 559 const PetscReal *varray; 560 const char *vecname; 561 PetscInt xs, xe, bs; 562 PetscBool useNatural; 563 int offset; 564 IS compIS; 565 PetscInt *csSize, *csID; 566 PetscInt numCS, set, csxs = 0; 567 PetscErrorCode ierr; 568 569 PetscFunctionBegin; 570 ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr); 571 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 572 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 573 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 574 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 575 if (useNatural) { 576 ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 577 ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 578 ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 579 } else { 580 vNatural = v; 581 } 582 /* Get the location of the variable in the exodus file */ 583 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 584 ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 585 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 586 /* Write local chunk of the result in the exodus file 587 exodus stores each component of a vector-valued field as a separate variable. 588 We assume that they are stored sequentially 589 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 590 but once the vector has been reordered to natural size, we cannot use the label informations 591 to figure out what to save where. */ 592 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 593 ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 594 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 595 for (set = 0; set < numCS; ++set) { 596 ex_block block; 597 598 block.id = csID[set]; 599 block.type = EX_ELEM_BLOCK; 600 PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 601 csSize[set] = block.num_entry; 602 } 603 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 604 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 605 if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 606 for (set = 0; set < numCS; set++) { 607 PetscInt csLocalSize, c; 608 609 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 610 local slice of zonal values: xs/bs,xm/bs-1 611 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 612 csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 613 if (bs == 1) { 614 ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 615 PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 616 ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 617 } else { 618 for (c = 0; c < bs; ++c) { 619 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 620 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 621 ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 622 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)])); 623 ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 624 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 625 } 626 } 627 csxs += csSize[set]; 628 } 629 ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 630 if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 631 if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 632 PetscFunctionReturn(0); 633 } 634 635 /* 636 VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 637 638 Collective on v 639 640 Input Parameters: 641 + v - The vector to be written 642 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 643 - step - the time step to read at (exodus steps are numbered starting from 1) 644 645 Notes: 646 The exodus result zonal variable index is obtained by comparing the Vec name and the 647 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 648 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 649 amongst all zonal variables. 650 651 Level: beginner 652 653 .keywords: mesh,ExodusII 654 .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal() 655 */ 656 PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 657 { 658 MPI_Comm comm; 659 PetscMPIInt size; 660 DM dm; 661 Vec vNatural, vComp; 662 PetscReal *varray; 663 const char *vecname; 664 PetscInt xs, xe, bs; 665 PetscBool useNatural; 666 int offset; 667 IS compIS; 668 PetscInt *csSize, *csID; 669 PetscInt numCS, set, csxs = 0; 670 PetscErrorCode ierr; 671 672 PetscFunctionBegin; 673 ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 674 ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 675 ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 676 ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 677 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 678 if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 679 else {vNatural = v;} 680 /* Get the location of the variable in the exodus file */ 681 ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 682 ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 683 if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 684 /* Read local chunk of the result in the exodus file 685 exodus stores each component of a vector-valued field as a separate variable. 686 We assume that they are stored sequentially 687 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 688 but once the vector has been reordered to natural size, we cannot use the label informations 689 to figure out what to save where. */ 690 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 691 ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 692 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 693 for (set = 0; set < numCS; ++set) { 694 ex_block block; 695 696 block.id = csID[set]; 697 block.type = EX_ELEM_BLOCK; 698 PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 699 csSize[set] = block.num_entry; 700 } 701 ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 702 ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 703 if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 704 for (set = 0; set < numCS; ++set) { 705 PetscInt csLocalSize, c; 706 707 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 708 local slice of zonal values: xs/bs,xm/bs-1 709 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 710 csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 711 if (bs == 1) { 712 ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 713 PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 714 ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 715 } else { 716 for (c = 0; c < bs; ++c) { 717 ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 718 ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 719 ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 720 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)])); 721 ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 722 ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 723 } 724 } 725 csxs += csSize[set]; 726 } 727 ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 728 if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 729 if (useNatural) { 730 ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 731 ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 732 ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 733 } 734 PetscFunctionReturn(0); 735 } 736 #endif 737 738 /*@C 739 DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 740 741 Collective on comm 742 743 Input Parameters: 744 + comm - The MPI communicator 745 . filename - The name of the ExodusII file 746 - interpolate - Create faces and edges in the mesh 747 748 Output Parameter: 749 . dm - The DM object representing the mesh 750 751 Level: beginner 752 753 .keywords: mesh,ExodusII 754 .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 755 @*/ 756 PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 757 { 758 PetscMPIInt rank; 759 PetscErrorCode ierr; 760 #if defined(PETSC_HAVE_EXODUSII) 761 int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 762 float version; 763 #endif 764 765 PetscFunctionBegin; 766 PetscValidCharPointer(filename, 2); 767 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 768 #if defined(PETSC_HAVE_EXODUSII) 769 if (!rank) { 770 exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 771 if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 772 } 773 ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 774 if (!rank) {PetscStackCallStandard(ex_close,(exoid));} 775 #else 776 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 777 #endif 778 PetscFunctionReturn(0); 779 } 780 781 /*@ 782 DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 783 784 Collective on comm 785 786 Input Parameters: 787 + comm - The MPI communicator 788 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 789 - interpolate - Create faces and edges in the mesh 790 791 Output Parameter: 792 . dm - The DM object representing the mesh 793 794 Level: beginner 795 796 .keywords: mesh,ExodusII 797 .seealso: DMPLEX, DMCreate() 798 @*/ 799 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 800 { 801 #if defined(PETSC_HAVE_EXODUSII) 802 PetscMPIInt num_proc, rank; 803 PetscSection coordSection; 804 Vec coordinates; 805 PetscScalar *coords; 806 PetscInt coordSize, v; 807 PetscErrorCode ierr; 808 /* Read from ex_get_init() */ 809 char title[PETSC_MAX_PATH_LEN+1]; 810 int dim = 0, numVertices = 0, numCells = 0; 811 int num_cs = 0, num_vs = 0, num_fs = 0; 812 #endif 813 814 PetscFunctionBegin; 815 #if defined(PETSC_HAVE_EXODUSII) 816 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 817 ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 818 ierr = DMCreate(comm, dm);CHKERRQ(ierr); 819 ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 820 /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 821 if (!rank) { 822 ierr = PetscMemzero(title,(PETSC_MAX_PATH_LEN+1)*sizeof(char));CHKERRQ(ierr); 823 PetscStackCallStandard(ex_get_init,(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs)); 824 if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 825 } 826 ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 827 ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 828 ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 829 ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 830 ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 831 832 /* Read cell sets information */ 833 if (!rank) { 834 PetscInt *cone; 835 int c, cs, c_loc, v, v_loc; 836 /* Read from ex_get_elem_blk_ids() */ 837 int *cs_id; 838 /* Read from ex_get_elem_block() */ 839 char buffer[PETSC_MAX_PATH_LEN+1]; 840 int num_cell_in_set, num_vertex_per_cell, num_attr; 841 /* Read from ex_get_elem_conn() */ 842 int *cs_connect; 843 844 /* Get cell sets IDs */ 845 ierr = PetscMalloc1(num_cs, &cs_id);CHKERRQ(ierr); 846 PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id)); 847 /* Read the cell set connectivity table and build mesh topology 848 EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 849 /* First set sizes */ 850 for (cs = 0, c = 0; cs < num_cs; cs++) { 851 PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr)); 852 for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 853 ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 854 } 855 } 856 ierr = DMSetUp(*dm);CHKERRQ(ierr); 857 for (cs = 0, c = 0; cs < num_cs; cs++) { 858 PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr)); 859 ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 860 PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL)); 861 /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 862 for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 863 for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 864 cone[v_loc] = cs_connect[v]+numCells-1; 865 } 866 if (dim == 3) { 867 /* Tetrahedra are inverted */ 868 if (num_vertex_per_cell == 4) { 869 PetscInt tmp = cone[0]; 870 cone[0] = cone[1]; 871 cone[1] = tmp; 872 } 873 /* Hexahedra are inverted */ 874 if (num_vertex_per_cell == 8) { 875 PetscInt tmp = cone[1]; 876 cone[1] = cone[3]; 877 cone[3] = tmp; 878 } 879 } 880 ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 881 ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 882 } 883 ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 884 } 885 ierr = PetscFree(cs_id);CHKERRQ(ierr); 886 } 887 ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 888 ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 889 if (interpolate) { 890 DM idm; 891 892 ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 893 ierr = DMDestroy(dm);CHKERRQ(ierr); 894 *dm = idm; 895 } 896 897 /* Create vertex set label */ 898 if (!rank && (num_vs > 0)) { 899 int vs, v; 900 /* Read from ex_get_node_set_ids() */ 901 int *vs_id; 902 /* Read from ex_get_node_set_param() */ 903 int num_vertex_in_set; 904 /* Read from ex_get_node_set() */ 905 int *vs_vertex_list; 906 907 /* Get vertex set ids */ 908 ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 909 PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id)); 910 for (vs = 0; vs < num_vs; ++vs) { 911 PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL)); 912 ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 913 PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL)); 914 for (v = 0; v < num_vertex_in_set; ++v) { 915 ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 916 } 917 ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 918 } 919 ierr = PetscFree(vs_id);CHKERRQ(ierr); 920 } 921 /* Read coordinates */ 922 ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 923 ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 924 ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 925 ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 926 for (v = numCells; v < numCells+numVertices; ++v) { 927 ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 928 ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 929 } 930 ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 931 ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 932 ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 933 ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 934 ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 935 ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 936 ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 937 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 938 if (!rank) { 939 PetscReal *x, *y, *z; 940 941 ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 942 PetscStackCallStandard(ex_get_coord,(exoid, x, y, z)); 943 if (dim > 0) { 944 for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 945 } 946 if (dim > 1) { 947 for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 948 } 949 if (dim > 2) { 950 for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 951 } 952 ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 953 } 954 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 955 ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 956 ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 957 958 /* Create side set label */ 959 if (!rank && interpolate && (num_fs > 0)) { 960 int fs, f, voff; 961 /* Read from ex_get_side_set_ids() */ 962 int *fs_id; 963 /* Read from ex_get_side_set_param() */ 964 int num_side_in_set; 965 /* Read from ex_get_side_set_node_list() */ 966 int *fs_vertex_count_list, *fs_vertex_list; 967 /* Read side set labels */ 968 char fs_name[MAX_STR_LENGTH+1]; 969 970 /* Get side set ids */ 971 ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 972 PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id)); 973 for (fs = 0; fs < num_fs; ++fs) { 974 PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL)); 975 ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 976 PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list)); 977 /* Get the specific name associated with this side set ID. */ 978 int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 979 for (f = 0, voff = 0; f < num_side_in_set; ++f) { 980 const PetscInt *faces = NULL; 981 PetscInt faceSize = fs_vertex_count_list[f], numFaces; 982 PetscInt faceVertices[4], v; 983 984 if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 985 for (v = 0; v < faceSize; ++v, ++voff) { 986 faceVertices[v] = fs_vertex_list[voff]+numCells-1; 987 } 988 ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 989 if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 990 ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 991 /* Only add the label if one has been detected for this side set. */ 992 if (!fs_name_err) { 993 ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 994 } 995 ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 996 } 997 ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 998 } 999 ierr = PetscFree(fs_id);CHKERRQ(ierr); 1000 } 1001 #else 1002 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1003 #endif 1004 PetscFunctionReturn(0); 1005 } 1006