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