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 #include <petsc/private/viewerimpl.h> 10 #include <petsc/private/viewerexodusiiimpl.h> 11 #if defined(PETSC_HAVE_EXODUSII) 12 /* 13 PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator. 14 15 Collective 16 17 Input Parameter: 18 . comm - the MPI communicator to share the ExodusII PetscViewer 19 20 Level: intermediate 21 22 Notes: 23 misses Fortran bindings 24 25 Notes: 26 Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return 27 an error code. The GLVIS PetscViewer is usually used in the form 28 $ XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm)); 29 30 .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()` 31 */ 32 PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm) 33 { 34 PetscViewer viewer; 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer); 39 if (ierr) { 40 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 41 PetscFunctionReturn(NULL); 42 } 43 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 44 if (ierr) { 45 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 46 PetscFunctionReturn(NULL); 47 } 48 PetscFunctionReturn(viewer); 49 } 50 51 static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer) 52 { 53 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data; 54 55 PetscFunctionBegin; 56 if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename)); 57 if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid: %d\n", exo->exoid)); 58 if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode: %d\n", exo->btype)); 59 if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order: %" PetscInt_FMT "\n", exo->order)); 60 PetscFunctionReturn(0); 61 } 62 63 static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscViewer v, PetscOptionItems *PetscOptionsObject) 64 { 65 PetscFunctionBegin; 66 PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options"); 67 PetscOptionsHeadEnd(); 68 PetscFunctionReturn(0); 69 } 70 71 static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer) 72 { 73 PetscFunctionBegin; 74 PetscFunctionReturn(0); 75 } 76 77 static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer) 78 { 79 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 80 81 PetscFunctionBegin; 82 if (exo->exoid >= 0) PetscCallExternal(ex_close, exo->exoid); 83 PetscCall(PetscFree(exo->filename)); 84 PetscCall(PetscFree(exo)); 85 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 86 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 87 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 88 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 89 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetId_C", NULL)); 90 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetOrder_C", NULL)); 91 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerSetOrder_C", NULL)); 92 PetscFunctionReturn(0); 93 } 94 95 static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[]) 96 { 97 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 98 PetscMPIInt rank; 99 int CPU_word_size, IO_word_size, EXO_mode; 100 MPI_Info mpi_info = MPI_INFO_NULL; 101 float EXO_version; 102 103 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 104 CPU_word_size = sizeof(PetscReal); 105 IO_word_size = sizeof(PetscReal); 106 107 PetscFunctionBegin; 108 if (exo->exoid >= 0) { 109 PetscCallExternal(ex_close, exo->exoid); 110 exo->exoid = -1; 111 } 112 if (exo->filename) PetscCall(PetscFree(exo->filename)); 113 PetscCall(PetscStrallocpy(name, &exo->filename)); 114 switch (exo->btype) { 115 case FILE_MODE_READ: 116 EXO_mode = EX_READ; 117 break; 118 case FILE_MODE_APPEND: 119 case FILE_MODE_UPDATE: 120 case FILE_MODE_APPEND_UPDATE: 121 /* Will fail if the file does not already exist */ 122 EXO_mode = EX_WRITE; 123 break; 124 case FILE_MODE_WRITE: 125 /* 126 exodus only allows writing geometry upon file creation, so we will let DMView create the file. 127 */ 128 PetscFunctionReturn(0); 129 break; 130 default: 131 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 132 } 133 #if defined(PETSC_USE_64BIT_INDICES) 134 EXO_mode += EX_ALL_INT64_API; 135 #endif 136 exo->exoid = ex_open_par(name, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PETSC_COMM_WORLD, mpi_info); 137 PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name); 138 PetscFunctionReturn(0); 139 } 140 141 static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name) 142 { 143 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 144 145 PetscFunctionBegin; 146 *name = exo->filename; 147 PetscFunctionReturn(0); 148 } 149 150 static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type) 151 { 152 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 153 154 PetscFunctionBegin; 155 exo->btype = type; 156 PetscFunctionReturn(0); 157 } 158 159 static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type) 160 { 161 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 162 163 PetscFunctionBegin; 164 *type = exo->btype; 165 PetscFunctionReturn(0); 166 } 167 168 static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid) 169 { 170 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 171 172 PetscFunctionBegin; 173 *exoid = exo->exoid; 174 PetscFunctionReturn(0); 175 } 176 177 static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order) 178 { 179 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 180 181 PetscFunctionBegin; 182 *order = exo->order; 183 PetscFunctionReturn(0); 184 } 185 186 static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order) 187 { 188 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 189 190 PetscFunctionBegin; 191 exo->order = order; 192 PetscFunctionReturn(0); 193 } 194 195 /*MC 196 PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file 197 198 .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`, 199 `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 200 201 Level: beginner 202 M*/ 203 204 PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v) 205 { 206 PetscViewer_ExodusII *exo; 207 208 PetscFunctionBegin; 209 PetscCall(PetscNew(&exo)); 210 211 v->data = (void *)exo; 212 v->ops->destroy = PetscViewerDestroy_ExodusII; 213 v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII; 214 v->ops->setup = PetscViewerSetUp_ExodusII; 215 v->ops->view = PetscViewerView_ExodusII; 216 v->ops->flush = 0; 217 exo->btype = (PetscFileMode)-1; 218 exo->filename = 0; 219 exo->exoid = -1; 220 221 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_ExodusII)); 222 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_ExodusII)); 223 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_ExodusII)); 224 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_ExodusII)); 225 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetId_C", PetscViewerExodusIIGetId_ExodusII)); 226 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerSetOrder_C", PetscViewerExodusIISetOrder_ExodusII)); 227 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetOrder_C", PetscViewerExodusIIGetOrder_ExodusII)); 228 PetscFunctionReturn(0); 229 } 230 231 /* 232 EXOGetVarIndex - Locate a result in an exodus file based on its name 233 234 Collective 235 236 Input Parameters: 237 + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 238 . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 239 - name - the name of the result 240 241 Output Parameters: 242 . varIndex - the location in the exodus file of the result 243 244 Notes: 245 The exodus variable index is obtained by comparing name and the 246 names of zonal variables declared in the exodus file. For instance if name is "V" 247 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 248 amongst all variables of type obj_type. 249 250 Level: beginner 251 252 .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 253 */ 254 PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 255 { 256 int num_vars, i, j; 257 char ext_name[MAX_STR_LENGTH + 1], var_name[MAX_STR_LENGTH + 1]; 258 const int num_suffix = 5; 259 char *suffix[5]; 260 PetscBool flg; 261 262 PetscFunctionBegin; 263 suffix[0] = (char *)""; 264 suffix[1] = (char *)"_X"; 265 suffix[2] = (char *)"_XX"; 266 suffix[3] = (char *)"_1"; 267 suffix[4] = (char *)"_11"; 268 269 *varIndex = -1; 270 PetscCallExternal(ex_get_variable_param, exoid, obj_type, &num_vars); 271 for (i = 0; i < num_vars; ++i) { 272 PetscCallExternal(ex_get_variable_name, exoid, obj_type, i + 1, var_name); 273 for (j = 0; j < num_suffix; ++j) { 274 PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH)); 275 PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH)); 276 PetscCall(PetscStrcasecmp(ext_name, var_name, &flg)); 277 if (flg) *varIndex = i + 1; 278 } 279 } 280 PetscFunctionReturn(0); 281 } 282 283 /* 284 DMView_PlexExodusII - Write a DM to disk in exodus format 285 286 Collective on dm 287 288 Input Parameters: 289 + dm - The dm to be written 290 . viewer - an exodusII viewer 291 292 Notes: 293 Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 294 consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 295 296 If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices) 297 will be written. 298 299 DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 300 on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 301 The order of the mesh shall be set using PetscViewerExodusIISetOrder 302 It should be extended to use PetscFE objects. 303 304 This function will only handle TRI, TET, QUAD, and HEX cells. 305 Level: beginner 306 307 .seealso: 308 */ 309 PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer) 310 { 311 enum ElemType { 312 TRI, 313 QUAD, 314 TET, 315 HEX 316 }; 317 MPI_Comm comm; 318 PetscInt degree; /* the order of the mesh */ 319 /* Connectivity Variables */ 320 PetscInt cellsNotInConnectivity; 321 /* Cell Sets */ 322 DMLabel csLabel; 323 IS csIS; 324 const PetscInt *csIdx; 325 PetscInt num_cs, cs; 326 enum ElemType *type; 327 PetscBool hasLabel; 328 /* Coordinate Variables */ 329 DM cdm; 330 PetscSection coordSection; 331 Vec coord; 332 PetscInt **nodes; 333 PetscInt depth, d, dim, skipCells = 0; 334 PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 335 PetscInt num_vs, num_fs; 336 PetscMPIInt rank, size; 337 const char *dmName; 338 PetscInt nodesTriP1[4] = {3, 0, 0, 0}; 339 PetscInt nodesTriP2[4] = {3, 3, 0, 0}; 340 PetscInt nodesQuadP1[4] = {4, 0, 0, 0}; 341 PetscInt nodesQuadP2[4] = {4, 4, 0, 1}; 342 PetscInt nodesTetP1[4] = {4, 0, 0, 0}; 343 PetscInt nodesTetP2[4] = {4, 6, 0, 0}; 344 PetscInt nodesHexP1[4] = {8, 0, 0, 0}; 345 PetscInt nodesHexP2[4] = {8, 12, 6, 1}; 346 int CPU_word_size, IO_word_size, EXO_mode; 347 float EXO_version; 348 349 PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 350 351 PetscFunctionBegin; 352 PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 353 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 354 PetscCallMPI(MPI_Comm_size(comm, &size)); 355 356 /* 357 Creating coordSection is a collective operation so we do it somewhat out of sequence 358 */ 359 PetscCall(PetscSectionCreate(comm, &coordSection)); 360 PetscCall(DMGetCoordinatesLocalSetUp(dm)); 361 if (rank == 0) { 362 switch (exo->btype) { 363 case FILE_MODE_READ: 364 case FILE_MODE_APPEND: 365 case FILE_MODE_UPDATE: 366 case FILE_MODE_APPEND_UPDATE: 367 /* exodusII does not allow writing geometry to an existing file */ 368 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename); 369 case FILE_MODE_WRITE: 370 /* Create an empty file if one already exists*/ 371 EXO_mode = EX_CLOBBER; 372 #if defined(PETSC_USE_64BIT_INDICES) 373 EXO_mode += EX_ALL_INT64_API; 374 #endif 375 CPU_word_size = sizeof(PetscReal); 376 IO_word_size = sizeof(PetscReal); 377 exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size); 378 PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename); 379 380 break; 381 default: 382 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 383 } 384 385 /* --- Get DM info --- */ 386 PetscCall(PetscObjectGetName((PetscObject)dm, &dmName)); 387 PetscCall(DMPlexGetDepth(dm, &depth)); 388 PetscCall(DMGetDimension(dm, &dim)); 389 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 390 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 391 PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 392 PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 393 PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 394 numCells = cEnd - cStart; 395 numEdges = eEnd - eStart; 396 numVertices = vEnd - vStart; 397 if (depth == 3) { 398 numFaces = fEnd - fStart; 399 } else { 400 numFaces = 0; 401 } 402 PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs)); 403 PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs)); 404 PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs)); 405 PetscCall(DMGetCoordinatesLocal(dm, &coord)); 406 PetscCall(DMGetCoordinateDM(dm, &cdm)); 407 if (num_cs > 0) { 408 PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel)); 409 PetscCall(DMLabelGetValueIS(csLabel, &csIS)); 410 PetscCall(ISGetIndices(csIS, &csIdx)); 411 } 412 PetscCall(PetscMalloc1(num_cs, &nodes)); 413 /* Set element type for each block and compute total number of nodes */ 414 PetscCall(PetscMalloc1(num_cs, &type)); 415 numNodes = numVertices; 416 417 PetscCall(PetscViewerExodusIIGetOrder(viewer, °ree)); 418 if (degree == 2) numNodes += numEdges; 419 cellsNotInConnectivity = numCells; 420 for (cs = 0; cs < num_cs; ++cs) { 421 IS stratumIS; 422 const PetscInt *cells; 423 PetscScalar *xyz = NULL; 424 PetscInt csSize, closureSize; 425 426 PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 427 PetscCall(ISGetIndices(stratumIS, &cells)); 428 PetscCall(ISGetSize(stratumIS, &csSize)); 429 PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 430 switch (dim) { 431 case 2: 432 if (closureSize == 3 * dim) { 433 type[cs] = TRI; 434 } else if (closureSize == 4 * dim) { 435 type[cs] = QUAD; 436 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim); 437 break; 438 case 3: 439 if (closureSize == 4 * dim) { 440 type[cs] = TET; 441 } else if (closureSize == 8 * dim) { 442 type[cs] = HEX; 443 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim); 444 break; 445 default: 446 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim); 447 } 448 if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize; 449 if ((degree == 2) && (type[cs] == HEX)) { 450 numNodes += csSize; 451 numNodes += numFaces; 452 } 453 PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 454 /* Set nodes and Element type */ 455 if (type[cs] == TRI) { 456 if (degree == 1) nodes[cs] = nodesTriP1; 457 else if (degree == 2) nodes[cs] = nodesTriP2; 458 } else if (type[cs] == QUAD) { 459 if (degree == 1) nodes[cs] = nodesQuadP1; 460 else if (degree == 2) nodes[cs] = nodesQuadP2; 461 } else if (type[cs] == TET) { 462 if (degree == 1) nodes[cs] = nodesTetP1; 463 else if (degree == 2) nodes[cs] = nodesTetP2; 464 } else if (type[cs] == HEX) { 465 if (degree == 1) nodes[cs] = nodesHexP1; 466 else if (degree == 2) nodes[cs] = nodesHexP2; 467 } 468 /* Compute the number of cells not in the connectivity table */ 469 cellsNotInConnectivity -= nodes[cs][3] * csSize; 470 471 PetscCall(ISRestoreIndices(stratumIS, &cells)); 472 PetscCall(ISDestroy(&stratumIS)); 473 } 474 if (num_cs > 0) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs); 475 /* --- Connectivity --- */ 476 for (cs = 0; cs < num_cs; ++cs) { 477 IS stratumIS; 478 const PetscInt *cells; 479 PetscInt *connect, off = 0; 480 PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 481 PetscInt csSize, c, connectSize, closureSize; 482 char *elem_type = NULL; 483 char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 484 char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 485 char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 486 char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 487 488 PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 489 PetscCall(ISGetIndices(stratumIS, &cells)); 490 PetscCall(ISGetSize(stratumIS, &csSize)); 491 /* Set Element type */ 492 if (type[cs] == TRI) { 493 if (degree == 1) elem_type = elem_type_tri3; 494 else if (degree == 2) elem_type = elem_type_tri6; 495 } else if (type[cs] == QUAD) { 496 if (degree == 1) elem_type = elem_type_quad4; 497 else if (degree == 2) elem_type = elem_type_quad9; 498 } else if (type[cs] == TET) { 499 if (degree == 1) elem_type = elem_type_tet4; 500 else if (degree == 2) elem_type = elem_type_tet10; 501 } else if (type[cs] == HEX) { 502 if (degree == 1) elem_type = elem_type_hex8; 503 else if (degree == 2) elem_type = elem_type_hex27; 504 } 505 connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 506 PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect)); 507 PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1); 508 /* Find number of vertices, edges, and faces in the closure */ 509 verticesInClosure = nodes[cs][0]; 510 if (depth > 1) { 511 if (dim == 2) { 512 PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure)); 513 } else if (dim == 3) { 514 PetscInt *closure = NULL; 515 516 PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure)); 517 PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 518 edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 519 PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 520 } 521 } 522 /* Get connectivity for each cell */ 523 for (c = 0; c < csSize; ++c) { 524 PetscInt *closure = NULL; 525 PetscInt temp, i; 526 527 PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 528 for (i = 0; i < connectSize; ++i) { 529 if (i < nodes[cs][0]) { /* Vertices */ 530 connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1; 531 connect[i + off] -= cellsNotInConnectivity; 532 } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */ 533 connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1; 534 if (nodes[cs][2] == 0) connect[i + off] -= numFaces; 535 connect[i + off] -= cellsNotInConnectivity; 536 } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */ 537 connect[i + off] = closure[0] + 1; 538 connect[i + off] -= skipCells; 539 } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */ 540 connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1; 541 connect[i + off] -= cellsNotInConnectivity; 542 } else { 543 connect[i + off] = -1; 544 } 545 } 546 /* Tetrahedra are inverted */ 547 if (type[cs] == TET) { 548 temp = connect[0 + off]; 549 connect[0 + off] = connect[1 + off]; 550 connect[1 + off] = temp; 551 if (degree == 2) { 552 temp = connect[5 + off]; 553 connect[5 + off] = connect[6 + off]; 554 connect[6 + off] = temp; 555 temp = connect[7 + off]; 556 connect[7 + off] = connect[8 + off]; 557 connect[8 + off] = temp; 558 } 559 } 560 /* Hexahedra are inverted */ 561 if (type[cs] == HEX) { 562 temp = connect[1 + off]; 563 connect[1 + off] = connect[3 + off]; 564 connect[3 + off] = temp; 565 if (degree == 2) { 566 temp = connect[8 + off]; 567 connect[8 + off] = connect[11 + off]; 568 connect[11 + off] = temp; 569 temp = connect[9 + off]; 570 connect[9 + off] = connect[10 + off]; 571 connect[10 + off] = temp; 572 temp = connect[16 + off]; 573 connect[16 + off] = connect[17 + off]; 574 connect[17 + off] = temp; 575 temp = connect[18 + off]; 576 connect[18 + off] = connect[19 + off]; 577 connect[19 + off] = temp; 578 579 temp = connect[12 + off]; 580 connect[12 + off] = connect[16 + off]; 581 connect[16 + off] = temp; 582 temp = connect[13 + off]; 583 connect[13 + off] = connect[17 + off]; 584 connect[17 + off] = temp; 585 temp = connect[14 + off]; 586 connect[14 + off] = connect[18 + off]; 587 connect[18 + off] = temp; 588 temp = connect[15 + off]; 589 connect[15 + off] = connect[19 + off]; 590 connect[19 + off] = temp; 591 592 temp = connect[23 + off]; 593 connect[23 + off] = connect[26 + off]; 594 connect[26 + off] = temp; 595 temp = connect[24 + off]; 596 connect[24 + off] = connect[25 + off]; 597 connect[25 + off] = temp; 598 temp = connect[25 + off]; 599 connect[25 + off] = connect[26 + off]; 600 connect[26 + off] = temp; 601 } 602 } 603 off += connectSize; 604 PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 605 } 606 PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0); 607 skipCells += (nodes[cs][3] == 0) * csSize; 608 PetscCall(PetscFree(connect)); 609 PetscCall(ISRestoreIndices(stratumIS, &cells)); 610 PetscCall(ISDestroy(&stratumIS)); 611 } 612 PetscCall(PetscFree(type)); 613 /* --- Coordinates --- */ 614 PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd)); 615 if (num_cs) { 616 for (d = 0; d < depth; ++d) { 617 PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 618 for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0)); 619 } 620 } 621 for (cs = 0; cs < num_cs; ++cs) { 622 IS stratumIS; 623 const PetscInt *cells; 624 PetscInt csSize, c; 625 626 PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 627 PetscCall(ISGetIndices(stratumIS, &cells)); 628 PetscCall(ISGetSize(stratumIS, &csSize)); 629 for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0)); 630 PetscCall(ISRestoreIndices(stratumIS, &cells)); 631 PetscCall(ISDestroy(&stratumIS)); 632 } 633 if (num_cs > 0) { 634 PetscCall(ISRestoreIndices(csIS, &csIdx)); 635 PetscCall(ISDestroy(&csIS)); 636 } 637 PetscCall(PetscFree(nodes)); 638 PetscCall(PetscSectionSetUp(coordSection)); 639 if (numNodes > 0) { 640 const char *coordNames[3] = {"x", "y", "z"}; 641 PetscScalar *closure, *cval; 642 PetscReal *coords; 643 PetscInt hasDof, n = 0; 644 645 /* There can't be more than 24 values in the closure of a point for the coord coordSection */ 646 PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure)); 647 PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord)); 648 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 649 for (p = pStart; p < pEnd; ++p) { 650 PetscCall(PetscSectionGetDof(coordSection, p, &hasDof)); 651 if (hasDof) { 652 PetscInt closureSize = 24, j; 653 654 PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure)); 655 for (d = 0; d < dim; ++d) { 656 cval[d] = 0.0; 657 for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d]; 658 coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize; 659 } 660 ++n; 661 } 662 } 663 PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]); 664 PetscCall(PetscFree3(coords, cval, closure)); 665 PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames); 666 } 667 668 /* --- Node Sets/Vertex Sets --- */ 669 DMHasLabel(dm, "Vertex Sets", &hasLabel); 670 if (hasLabel) { 671 PetscInt i, vs, vsSize; 672 const PetscInt *vsIdx, *vertices; 673 PetscInt *nodeList; 674 IS vsIS, stratumIS; 675 DMLabel vsLabel; 676 PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel)); 677 PetscCall(DMLabelGetValueIS(vsLabel, &vsIS)); 678 PetscCall(ISGetIndices(vsIS, &vsIdx)); 679 for (vs = 0; vs < num_vs; ++vs) { 680 PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS)); 681 PetscCall(ISGetIndices(stratumIS, &vertices)); 682 PetscCall(ISGetSize(stratumIS, &vsSize)); 683 PetscCall(PetscMalloc1(vsSize, &nodeList)); 684 for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1; 685 PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0); 686 PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL); 687 PetscCall(ISRestoreIndices(stratumIS, &vertices)); 688 PetscCall(ISDestroy(&stratumIS)); 689 PetscCall(PetscFree(nodeList)); 690 } 691 PetscCall(ISRestoreIndices(vsIS, &vsIdx)); 692 PetscCall(ISDestroy(&vsIS)); 693 } 694 /* --- Side Sets/Face Sets --- */ 695 PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel)); 696 if (hasLabel) { 697 PetscInt i, j, fs, fsSize; 698 const PetscInt *fsIdx, *faces; 699 IS fsIS, stratumIS; 700 DMLabel fsLabel; 701 PetscInt numPoints, *points; 702 PetscInt elem_list_size = 0; 703 PetscInt *elem_list, *elem_ind, *side_list; 704 705 PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel)); 706 /* Compute size of Node List and Element List */ 707 PetscCall(DMLabelGetValueIS(fsLabel, &fsIS)); 708 PetscCall(ISGetIndices(fsIS, &fsIdx)); 709 for (fs = 0; fs < num_fs; ++fs) { 710 PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 711 PetscCall(ISGetSize(stratumIS, &fsSize)); 712 elem_list_size += fsSize; 713 PetscCall(ISDestroy(&stratumIS)); 714 } 715 PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list)); 716 elem_ind[0] = 0; 717 for (fs = 0; fs < num_fs; ++fs) { 718 PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 719 PetscCall(ISGetIndices(stratumIS, &faces)); 720 PetscCall(ISGetSize(stratumIS, &fsSize)); 721 /* Set Parameters */ 722 PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0); 723 /* Indices */ 724 if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize; 725 726 for (i = 0; i < fsSize; ++i) { 727 /* Element List */ 728 points = NULL; 729 PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 730 elem_list[elem_ind[fs] + i] = points[2] + 1; 731 PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 732 733 /* Side List */ 734 points = NULL; 735 PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 736 for (j = 1; j < numPoints; ++j) { 737 if (points[j * 2] == faces[i]) break; 738 } 739 /* Convert HEX sides */ 740 if (numPoints == 27) { 741 if (j == 1) { 742 j = 5; 743 } else if (j == 2) { 744 j = 6; 745 } else if (j == 3) { 746 j = 1; 747 } else if (j == 4) { 748 j = 3; 749 } else if (j == 5) { 750 j = 2; 751 } else if (j == 6) { 752 j = 4; 753 } 754 } 755 /* Convert TET sides */ 756 if (numPoints == 15) { 757 --j; 758 if (j == 0) j = 4; 759 } 760 side_list[elem_ind[fs] + i] = j; 761 PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 762 } 763 PetscCall(ISRestoreIndices(stratumIS, &faces)); 764 PetscCall(ISDestroy(&stratumIS)); 765 } 766 PetscCall(ISRestoreIndices(fsIS, &fsIdx)); 767 PetscCall(ISDestroy(&fsIS)); 768 769 /* Put side sets */ 770 for (fs = 0; fs < num_fs; ++fs) PetscCallExternal(ex_put_set, exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]); 771 PetscCall(PetscFree3(elem_ind, elem_list, side_list)); 772 } 773 /* 774 close the exodus file 775 */ 776 ex_close(exo->exoid); 777 exo->exoid = -1; 778 } 779 PetscCall(PetscSectionDestroy(&coordSection)); 780 781 /* 782 reopen the file in parallel 783 */ 784 EXO_mode = EX_WRITE; 785 #if defined(PETSC_USE_64BIT_INDICES) 786 EXO_mode += EX_ALL_INT64_API; 787 #endif 788 CPU_word_size = sizeof(PetscReal); 789 IO_word_size = sizeof(PetscReal); 790 exo->exoid = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL); 791 PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename); 792 PetscFunctionReturn(0); 793 } 794 795 /* 796 VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 797 798 Collective on v 799 800 Input Parameters: 801 + v - The vector to be written 802 - viewer - The PetscViewerExodusII viewer associate to an exodus file 803 804 Notes: 805 The exodus result variable index is obtained by comparing the Vec name and the 806 names of variables declared in the exodus file. For instance for a Vec named "V" 807 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 808 amongst all variables. 809 In the event where a nodal and zonal variable both match, the function will return an error instead of 810 possibly corrupting the file 811 812 Level: beginner 813 814 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 815 @*/ 816 PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer) 817 { 818 DM dm; 819 MPI_Comm comm; 820 PetscMPIInt rank; 821 822 int exoid, offsetN = 0, offsetZ = 0; 823 const char *vecname; 824 PetscInt step; 825 826 PetscFunctionBegin; 827 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 828 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 829 PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 830 PetscCall(VecGetDM(v, &dm)); 831 PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 832 833 PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 834 PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 835 PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 836 PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 837 if (offsetN > 0) { 838 PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 839 } else if (offsetZ > 0) { 840 PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 841 } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 842 PetscFunctionReturn(0); 843 } 844 845 /* 846 VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 847 848 Collective on v 849 850 Input Parameters: 851 + v - The vector to be written 852 - viewer - The PetscViewerExodusII viewer associate to an exodus file 853 854 Notes: 855 The exodus result variable index is obtained by comparing the Vec name and the 856 names of variables declared in the exodus file. For instance for a Vec named "V" 857 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 858 amongst all variables. 859 In the event where a nodal and zonal variable both match, the function will return an error instead of 860 possibly corrupting the file 861 862 Level: beginner 863 864 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 865 @*/ 866 PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer) 867 { 868 DM dm; 869 MPI_Comm comm; 870 PetscMPIInt rank; 871 872 int exoid, offsetN = 0, offsetZ = 0; 873 const char *vecname; 874 PetscInt step; 875 876 PetscFunctionBegin; 877 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 878 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 879 PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 880 PetscCall(VecGetDM(v, &dm)); 881 PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 882 883 PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 884 PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 885 PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 886 PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 887 if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 888 else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 889 else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 890 PetscFunctionReturn(0); 891 } 892 893 /* 894 VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 895 896 Collective on v 897 898 Input Parameters: 899 + v - The vector to be written 900 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 901 . step - the time step to write at (exodus steps are numbered starting from 1) 902 - offset - the location of the variable in the file 903 904 Notes: 905 The exodus result nodal variable index is obtained by comparing the Vec name and the 906 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 907 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 908 amongst all nodal variables. 909 910 Level: beginner 911 912 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 913 @*/ 914 PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 915 { 916 MPI_Comm comm; 917 PetscMPIInt size; 918 DM dm; 919 Vec vNatural, vComp; 920 const PetscScalar *varray; 921 PetscInt xs, xe, bs; 922 PetscBool useNatural; 923 924 PetscFunctionBegin; 925 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 926 PetscCallMPI(MPI_Comm_size(comm, &size)); 927 PetscCall(VecGetDM(v, &dm)); 928 PetscCall(DMGetUseNatural(dm, &useNatural)); 929 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 930 if (useNatural) { 931 PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 932 PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 933 PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 934 } else { 935 vNatural = v; 936 } 937 938 /* Write local chunk of the result in the exodus file 939 exodus stores each component of a vector-valued field as a separate variable. 940 We assume that they are stored sequentially */ 941 PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 942 PetscCall(VecGetBlockSize(vNatural, &bs)); 943 if (bs == 1) { 944 PetscCall(VecGetArrayRead(vNatural, &varray)); 945 PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 946 PetscCall(VecRestoreArrayRead(vNatural, &varray)); 947 } else { 948 IS compIS; 949 PetscInt c; 950 951 PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 952 for (c = 0; c < bs; ++c) { 953 PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 954 PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 955 PetscCall(VecGetArrayRead(vComp, &varray)); 956 PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 957 PetscCall(VecRestoreArrayRead(vComp, &varray)); 958 PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 959 } 960 PetscCall(ISDestroy(&compIS)); 961 } 962 if (useNatural) PetscCall(VecDestroy(&vNatural)); 963 PetscFunctionReturn(0); 964 } 965 966 /* 967 VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 968 969 Collective on v 970 971 Input Parameters: 972 + v - The vector to be written 973 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 974 . step - the time step to read at (exodus steps are numbered starting from 1) 975 - offset - the location of the variable in the file 976 977 Notes: 978 The exodus result nodal variable index is obtained by comparing the Vec name and the 979 names of nodal variables declared in the exodus file. For instance for a Vec named "V" 980 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 981 amongst all nodal variables. 982 983 Level: beginner 984 985 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 986 */ 987 PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 988 { 989 MPI_Comm comm; 990 PetscMPIInt size; 991 DM dm; 992 Vec vNatural, vComp; 993 PetscScalar *varray; 994 PetscInt xs, xe, bs; 995 PetscBool useNatural; 996 997 PetscFunctionBegin; 998 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 999 PetscCallMPI(MPI_Comm_size(comm, &size)); 1000 PetscCall(VecGetDM(v, &dm)); 1001 PetscCall(DMGetUseNatural(dm, &useNatural)); 1002 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1003 if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1004 else vNatural = v; 1005 1006 /* Read local chunk from the file */ 1007 PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 1008 PetscCall(VecGetBlockSize(vNatural, &bs)); 1009 if (bs == 1) { 1010 PetscCall(VecGetArray(vNatural, &varray)); 1011 PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 1012 PetscCall(VecRestoreArray(vNatural, &varray)); 1013 } else { 1014 IS compIS; 1015 PetscInt c; 1016 1017 PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 1018 for (c = 0; c < bs; ++c) { 1019 PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 1020 PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 1021 PetscCall(VecGetArray(vComp, &varray)); 1022 PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 1023 PetscCall(VecRestoreArray(vComp, &varray)); 1024 PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 1025 } 1026 PetscCall(ISDestroy(&compIS)); 1027 } 1028 if (useNatural) { 1029 PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 1030 PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1031 PetscCall(VecDestroy(&vNatural)); 1032 } 1033 PetscFunctionReturn(0); 1034 } 1035 1036 /* 1037 VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 1038 1039 Collective on v 1040 1041 Input Parameters: 1042 + v - The vector to be written 1043 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 1044 . step - the time step to write at (exodus steps are numbered starting from 1) 1045 - offset - the location of the variable in the file 1046 1047 Notes: 1048 The exodus result zonal variable index is obtained by comparing the Vec name and the 1049 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 1050 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 1051 amongst all zonal variables. 1052 1053 Level: beginner 1054 1055 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 1056 */ 1057 PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1058 { 1059 MPI_Comm comm; 1060 PetscMPIInt size; 1061 DM dm; 1062 Vec vNatural, vComp; 1063 const PetscScalar *varray; 1064 PetscInt xs, xe, bs; 1065 PetscBool useNatural; 1066 IS compIS; 1067 PetscInt *csSize, *csID; 1068 PetscInt numCS, set, csxs = 0; 1069 1070 PetscFunctionBegin; 1071 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 1072 PetscCallMPI(MPI_Comm_size(comm, &size)); 1073 PetscCall(VecGetDM(v, &dm)); 1074 PetscCall(DMGetUseNatural(dm, &useNatural)); 1075 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1076 if (useNatural) { 1077 PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1078 PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 1079 PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 1080 } else { 1081 vNatural = v; 1082 } 1083 1084 /* Write local chunk of the result in the exodus file 1085 exodus stores each component of a vector-valued field as a separate variable. 1086 We assume that they are stored sequentially 1087 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1088 but once the vector has been reordered to natural size, we cannot use the label information 1089 to figure out what to save where. */ 1090 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 1091 PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1092 PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 1093 for (set = 0; set < numCS; ++set) { 1094 ex_block block; 1095 1096 block.id = csID[set]; 1097 block.type = EX_ELEM_BLOCK; 1098 PetscCallExternal(ex_get_block_param, exoid, &block); 1099 csSize[set] = block.num_entry; 1100 } 1101 PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 1102 PetscCall(VecGetBlockSize(vNatural, &bs)); 1103 if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 1104 for (set = 0; set < numCS; set++) { 1105 PetscInt csLocalSize, c; 1106 1107 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 1108 local slice of zonal values: xs/bs,xm/bs-1 1109 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 1110 csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 1111 if (bs == 1) { 1112 PetscCall(VecGetArrayRead(vNatural, &varray)); 1113 PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]); 1114 PetscCall(VecRestoreArrayRead(vNatural, &varray)); 1115 } else { 1116 for (c = 0; c < bs; ++c) { 1117 PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 1118 PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 1119 PetscCall(VecGetArrayRead(vComp, &varray)); 1120 PetscCallExternal(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)]); 1121 PetscCall(VecRestoreArrayRead(vComp, &varray)); 1122 PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 1123 } 1124 } 1125 csxs += csSize[set]; 1126 } 1127 PetscCall(PetscFree2(csID, csSize)); 1128 if (bs > 1) PetscCall(ISDestroy(&compIS)); 1129 if (useNatural) PetscCall(VecDestroy(&vNatural)); 1130 PetscFunctionReturn(0); 1131 } 1132 1133 /* 1134 VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 1135 1136 Collective on v 1137 1138 Input Parameters: 1139 + v - The vector to be written 1140 . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 1141 . step - the time step to read at (exodus steps are numbered starting from 1) 1142 - offset - the location of the variable in the file 1143 1144 Notes: 1145 The exodus result zonal variable index is obtained by comparing the Vec name and the 1146 names of zonal variables declared in the exodus file. For instance for a Vec named "V" 1147 the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 1148 amongst all zonal variables. 1149 1150 Level: beginner 1151 1152 .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 1153 */ 1154 PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1155 { 1156 MPI_Comm comm; 1157 PetscMPIInt size; 1158 DM dm; 1159 Vec vNatural, vComp; 1160 PetscScalar *varray; 1161 PetscInt xs, xe, bs; 1162 PetscBool useNatural; 1163 IS compIS; 1164 PetscInt *csSize, *csID; 1165 PetscInt numCS, set, csxs = 0; 1166 1167 PetscFunctionBegin; 1168 PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 1169 PetscCallMPI(MPI_Comm_size(comm, &size)); 1170 PetscCall(VecGetDM(v, &dm)); 1171 PetscCall(DMGetUseNatural(dm, &useNatural)); 1172 useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1173 if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1174 else vNatural = v; 1175 1176 /* Read local chunk of the result in the exodus file 1177 exodus stores each component of a vector-valued field as a separate variable. 1178 We assume that they are stored sequentially 1179 Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1180 but once the vector has been reordered to natural size, we cannot use the label information 1181 to figure out what to save where. */ 1182 numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 1183 PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1184 PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 1185 for (set = 0; set < numCS; ++set) { 1186 ex_block block; 1187 1188 block.id = csID[set]; 1189 block.type = EX_ELEM_BLOCK; 1190 PetscCallExternal(ex_get_block_param, exoid, &block); 1191 csSize[set] = block.num_entry; 1192 } 1193 PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 1194 PetscCall(VecGetBlockSize(vNatural, &bs)); 1195 if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 1196 for (set = 0; set < numCS; ++set) { 1197 PetscInt csLocalSize, c; 1198 1199 /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 1200 local slice of zonal values: xs/bs,xm/bs-1 1201 intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 1202 csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 1203 if (bs == 1) { 1204 PetscCall(VecGetArray(vNatural, &varray)); 1205 PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]); 1206 PetscCall(VecRestoreArray(vNatural, &varray)); 1207 } else { 1208 for (c = 0; c < bs; ++c) { 1209 PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 1210 PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 1211 PetscCall(VecGetArray(vComp, &varray)); 1212 PetscCallExternal(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)]); 1213 PetscCall(VecRestoreArray(vComp, &varray)); 1214 PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 1215 } 1216 } 1217 csxs += csSize[set]; 1218 } 1219 PetscCall(PetscFree2(csID, csSize)); 1220 if (bs > 1) PetscCall(ISDestroy(&compIS)); 1221 if (useNatural) { 1222 PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 1223 PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1224 PetscCall(VecDestroy(&vNatural)); 1225 } 1226 PetscFunctionReturn(0); 1227 } 1228 #endif 1229 1230 /*@ 1231 PetscViewerExodusIIGetId - Get the file id of the ExodusII file 1232 1233 Logically Collective on PetscViewer 1234 1235 Input Parameter: 1236 . viewer - the PetscViewer 1237 1238 Output Parameter: 1239 . exoid - The ExodusII file id 1240 1241 Level: intermediate 1242 1243 .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 1244 @*/ 1245 PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) 1246 { 1247 PetscFunctionBegin; 1248 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1249 PetscTryMethod(viewer, "PetscViewerGetId_C", (PetscViewer, int *), (viewer, exoid)); 1250 PetscFunctionReturn(0); 1251 } 1252 1253 /*@ 1254 PetscViewerExodusIISetOrder - Set the elements order in the exodusII file. 1255 1256 Collective 1257 1258 Input Parameters: 1259 + viewer - the viewer 1260 - order - elements order 1261 1262 Output Parameter: 1263 1264 Level: beginner 1265 1266 Note: 1267 1268 .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 1269 @*/ 1270 PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order) 1271 { 1272 PetscFunctionBegin; 1273 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1274 PetscTryMethod(viewer, "PetscViewerSetOrder_C", (PetscViewer, PetscInt), (viewer, order)); 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /*@ 1279 PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file. 1280 1281 Collective 1282 1283 Input Parameters: 1284 + viewer - the viewer 1285 - order - elements order 1286 1287 Output Parameter: 1288 1289 Level: beginner 1290 1291 Note: 1292 1293 .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 1294 @*/ 1295 PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order) 1296 { 1297 PetscFunctionBegin; 1298 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1299 PetscTryMethod(viewer, "PetscViewerGetOrder_C", (PetscViewer, PetscInt *), (viewer, order)); 1300 PetscFunctionReturn(0); 1301 } 1302 1303 /*@C 1304 PetscViewerExodusIIOpen - Opens a file for ExodusII input/output. 1305 1306 Collective 1307 1308 Input Parameters: 1309 + comm - MPI communicator 1310 . name - name of file 1311 - type - type of file 1312 $ FILE_MODE_WRITE - create new file for binary output 1313 $ FILE_MODE_READ - open existing file for binary input 1314 $ FILE_MODE_APPEND - open existing file for binary output 1315 1316 Output Parameter: 1317 . exo - PetscViewer for Exodus II input/output to use with the specified file 1318 1319 Level: beginner 1320 1321 Note: 1322 This PetscViewer should be destroyed with PetscViewerDestroy(). 1323 1324 .seealso: `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 1325 `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 1326 @*/ 1327 PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) 1328 { 1329 PetscFunctionBegin; 1330 PetscCall(PetscViewerCreate(comm, exo)); 1331 PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII)); 1332 PetscCall(PetscViewerFileSetMode(*exo, type)); 1333 PetscCall(PetscViewerFileSetName(*exo, name)); 1334 PetscCall(PetscViewerSetFromOptions(*exo)); 1335 PetscFunctionReturn(0); 1336 } 1337 1338 /*@C 1339 DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 1340 1341 Collective 1342 1343 Input Parameters: 1344 + comm - The MPI communicator 1345 . filename - The name of the ExodusII file 1346 - interpolate - Create faces and edges in the mesh 1347 1348 Output Parameter: 1349 . dm - The DM object representing the mesh 1350 1351 Level: beginner 1352 1353 .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()` 1354 @*/ 1355 PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1356 { 1357 PetscMPIInt rank; 1358 #if defined(PETSC_HAVE_EXODUSII) 1359 int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 1360 float version; 1361 #endif 1362 1363 PetscFunctionBegin; 1364 PetscValidCharPointer(filename, 2); 1365 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1366 #if defined(PETSC_HAVE_EXODUSII) 1367 if (rank == 0) { 1368 exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 1369 PetscCheck(exoid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 1370 } 1371 PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm)); 1372 if (rank == 0) PetscCallExternal(ex_close, exoid); 1373 PetscFunctionReturn(0); 1374 #else 1375 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1376 #endif 1377 } 1378 1379 #if defined(PETSC_HAVE_EXODUSII) 1380 static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct) 1381 { 1382 PetscBool flg; 1383 1384 PetscFunctionBegin; 1385 *ct = DM_POLYTOPE_UNKNOWN; 1386 PetscCall(PetscStrcmp(elem_type, "TRI", &flg)); 1387 if (flg) { 1388 *ct = DM_POLYTOPE_TRIANGLE; 1389 goto done; 1390 } 1391 PetscCall(PetscStrcmp(elem_type, "TRI3", &flg)); 1392 if (flg) { 1393 *ct = DM_POLYTOPE_TRIANGLE; 1394 goto done; 1395 } 1396 PetscCall(PetscStrcmp(elem_type, "QUAD", &flg)); 1397 if (flg) { 1398 *ct = DM_POLYTOPE_QUADRILATERAL; 1399 goto done; 1400 } 1401 PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg)); 1402 if (flg) { 1403 *ct = DM_POLYTOPE_QUADRILATERAL; 1404 goto done; 1405 } 1406 PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg)); 1407 if (flg) { 1408 *ct = DM_POLYTOPE_QUADRILATERAL; 1409 goto done; 1410 } 1411 PetscCall(PetscStrcmp(elem_type, "TETRA", &flg)); 1412 if (flg) { 1413 *ct = DM_POLYTOPE_TETRAHEDRON; 1414 goto done; 1415 } 1416 PetscCall(PetscStrcmp(elem_type, "TET4", &flg)); 1417 if (flg) { 1418 *ct = DM_POLYTOPE_TETRAHEDRON; 1419 goto done; 1420 } 1421 PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg)); 1422 if (flg) { 1423 *ct = DM_POLYTOPE_TRI_PRISM; 1424 goto done; 1425 } 1426 PetscCall(PetscStrcmp(elem_type, "HEX", &flg)); 1427 if (flg) { 1428 *ct = DM_POLYTOPE_HEXAHEDRON; 1429 goto done; 1430 } 1431 PetscCall(PetscStrcmp(elem_type, "HEX8", &flg)); 1432 if (flg) { 1433 *ct = DM_POLYTOPE_HEXAHEDRON; 1434 goto done; 1435 } 1436 PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg)); 1437 if (flg) { 1438 *ct = DM_POLYTOPE_HEXAHEDRON; 1439 goto done; 1440 } 1441 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type); 1442 done: 1443 PetscFunctionReturn(0); 1444 } 1445 #endif 1446 1447 /*@ 1448 DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 1449 1450 Collective 1451 1452 Input Parameters: 1453 + comm - The MPI communicator 1454 . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 1455 - interpolate - Create faces and edges in the mesh 1456 1457 Output Parameter: 1458 . dm - The DM object representing the mesh 1459 1460 Level: beginner 1461 1462 .seealso: `DMPLEX`, `DMCreate()` 1463 @*/ 1464 PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 1465 { 1466 #if defined(PETSC_HAVE_EXODUSII) 1467 PetscMPIInt num_proc, rank; 1468 DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL; 1469 PetscSection coordSection; 1470 Vec coordinates; 1471 PetscScalar *coords; 1472 PetscInt coordSize, v; 1473 /* Read from ex_get_init() */ 1474 char title[PETSC_MAX_PATH_LEN + 1]; 1475 int dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0; 1476 int num_cs = 0, num_vs = 0, num_fs = 0; 1477 #endif 1478 1479 PetscFunctionBegin; 1480 #if defined(PETSC_HAVE_EXODUSII) 1481 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1482 PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 1483 PetscCall(DMCreate(comm, dm)); 1484 PetscCall(DMSetType(*dm, DMPLEX)); 1485 /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */ 1486 if (rank == 0) { 1487 PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1)); 1488 PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 1489 PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set"); 1490 } 1491 PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm)); 1492 PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 1493 PetscCall(PetscObjectSetName((PetscObject)*dm, title)); 1494 PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 1495 /* We do not want this label automatically computed, instead we compute it here */ 1496 PetscCall(DMCreateLabel(*dm, "celltype")); 1497 1498 /* Read cell sets information */ 1499 if (rank == 0) { 1500 PetscInt *cone; 1501 int c, cs, ncs, c_loc, v, v_loc; 1502 /* Read from ex_get_elem_blk_ids() */ 1503 int *cs_id, *cs_order; 1504 /* Read from ex_get_elem_block() */ 1505 char buffer[PETSC_MAX_PATH_LEN + 1]; 1506 int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr; 1507 /* Read from ex_get_elem_conn() */ 1508 int *cs_connect; 1509 1510 /* Get cell sets IDs */ 1511 PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order)); 1512 PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id); 1513 /* Read the cell set connectivity table and build mesh topology 1514 EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 1515 /* Check for a hybrid mesh */ 1516 for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) { 1517 DMPolytopeType ct; 1518 char elem_type[PETSC_MAX_PATH_LEN]; 1519 1520 PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1521 PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 1522 PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 1523 dim = PetscMax(dim, DMPolytopeTypeGetDim(ct)); 1524 PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 1525 switch (ct) { 1526 case DM_POLYTOPE_TRI_PRISM: 1527 cs_order[cs] = cs; 1528 ++num_hybrid; 1529 break; 1530 default: 1531 for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1]; 1532 cs_order[cs - num_hybrid] = cs; 1533 } 1534 } 1535 /* First set sizes */ 1536 for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1537 DMPolytopeType ct; 1538 char elem_type[PETSC_MAX_PATH_LEN]; 1539 const PetscInt cs = cs_order[ncs]; 1540 1541 PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1542 PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 1543 PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 1544 PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 1545 for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1546 PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell)); 1547 PetscCall(DMPlexSetCellType(*dm, c, ct)); 1548 } 1549 } 1550 for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 1551 PetscCall(DMSetUp(*dm)); 1552 for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1553 const PetscInt cs = cs_order[ncs]; 1554 PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 1555 PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone)); 1556 PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL); 1557 /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 1558 for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1559 DMPolytopeType ct; 1560 1561 for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1; 1562 PetscCall(DMPlexGetCellType(*dm, c, &ct)); 1563 PetscCall(DMPlexInvertCell(ct, cone)); 1564 PetscCall(DMPlexSetCone(*dm, c, cone)); 1565 PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs])); 1566 } 1567 PetscCall(PetscFree2(cs_connect, cone)); 1568 } 1569 PetscCall(PetscFree2(cs_id, cs_order)); 1570 } 1571 { 1572 PetscInt ints[] = {dim, dimEmbed}; 1573 1574 PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm)); 1575 PetscCall(DMSetDimension(*dm, ints[0])); 1576 PetscCall(DMSetCoordinateDim(*dm, ints[1])); 1577 dim = ints[0]; 1578 dimEmbed = ints[1]; 1579 } 1580 PetscCall(DMPlexSymmetrize(*dm)); 1581 PetscCall(DMPlexStratify(*dm)); 1582 if (interpolate) { 1583 DM idm; 1584 1585 PetscCall(DMPlexInterpolate(*dm, &idm)); 1586 PetscCall(DMDestroy(dm)); 1587 *dm = idm; 1588 } 1589 1590 /* Create vertex set label */ 1591 if (rank == 0 && (num_vs > 0)) { 1592 int vs, v; 1593 /* Read from ex_get_node_set_ids() */ 1594 int *vs_id; 1595 /* Read from ex_get_node_set_param() */ 1596 int num_vertex_in_set; 1597 /* Read from ex_get_node_set() */ 1598 int *vs_vertex_list; 1599 1600 /* Get vertex set ids */ 1601 PetscCall(PetscMalloc1(num_vs, &vs_id)); 1602 PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id); 1603 for (vs = 0; vs < num_vs; ++vs) { 1604 PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 1605 PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list)); 1606 PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL); 1607 for (v = 0; v < num_vertex_in_set; ++v) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v] + numCells - 1, vs_id[vs])); 1608 PetscCall(PetscFree(vs_vertex_list)); 1609 } 1610 PetscCall(PetscFree(vs_id)); 1611 } 1612 /* Read coordinates */ 1613 PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 1614 PetscCall(PetscSectionSetNumFields(coordSection, 1)); 1615 PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed)); 1616 PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 1617 for (v = numCells; v < numCells + numVertices; ++v) { 1618 PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed)); 1619 PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed)); 1620 } 1621 PetscCall(PetscSectionSetUp(coordSection)); 1622 PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 1623 PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 1624 PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 1625 PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 1626 PetscCall(VecSetBlockSize(coordinates, dimEmbed)); 1627 PetscCall(VecSetType(coordinates, VECSTANDARD)); 1628 PetscCall(VecGetArray(coordinates, &coords)); 1629 if (rank == 0) { 1630 PetscReal *x, *y, *z; 1631 1632 PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z)); 1633 PetscCallExternal(ex_get_coord, exoid, x, y, z); 1634 if (dimEmbed > 0) { 1635 for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v]; 1636 } 1637 if (dimEmbed > 1) { 1638 for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v]; 1639 } 1640 if (dimEmbed > 2) { 1641 for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v]; 1642 } 1643 PetscCall(PetscFree3(x, y, z)); 1644 } 1645 PetscCall(VecRestoreArray(coordinates, &coords)); 1646 PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 1647 PetscCall(VecDestroy(&coordinates)); 1648 1649 /* Create side set label */ 1650 if (rank == 0 && interpolate && (num_fs > 0)) { 1651 int fs, f, voff; 1652 /* Read from ex_get_side_set_ids() */ 1653 int *fs_id; 1654 /* Read from ex_get_side_set_param() */ 1655 int num_side_in_set; 1656 /* Read from ex_get_side_set_node_list() */ 1657 int *fs_vertex_count_list, *fs_vertex_list; 1658 /* Read side set labels */ 1659 char fs_name[MAX_STR_LENGTH + 1]; 1660 size_t fs_name_len; 1661 1662 /* Get side set ids */ 1663 PetscCall(PetscMalloc1(num_fs, &fs_id)); 1664 PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id); 1665 for (fs = 0; fs < num_fs; ++fs) { 1666 PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL); 1667 PetscCall(PetscMalloc2(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list)); 1668 PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list); 1669 /* Get the specific name associated with this side set ID. */ 1670 int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1671 if (!fs_name_err) { 1672 PetscCall(PetscStrlen(fs_name, &fs_name_len)); 1673 if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH)); 1674 } 1675 for (f = 0, voff = 0; f < num_side_in_set; ++f) { 1676 const PetscInt *faces = NULL; 1677 PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1678 PetscInt faceVertices[4], v; 1679 1680 PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize); 1681 for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1; 1682 PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1683 PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces); 1684 PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs])); 1685 /* Only add the label if one has been detected for this side set. */ 1686 if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs])); 1687 PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1688 } 1689 PetscCall(PetscFree2(fs_vertex_count_list, fs_vertex_list)); 1690 } 1691 PetscCall(PetscFree(fs_id)); 1692 } 1693 1694 { /* Create Cell/Face/Vertex Sets labels at all processes */ 1695 enum { 1696 n = 3 1697 }; 1698 PetscBool flag[n]; 1699 1700 flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 1701 flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 1702 flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 1703 PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 1704 if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 1705 if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 1706 if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 1707 } 1708 PetscFunctionReturn(0); 1709 #else 1710 SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1711 #endif 1712 } 1713