1552f7358SJed Brown #define PETSCDM_DLL 2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 3552f7358SJed Brown 4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII) 5552f7358SJed Brown #include <netcdf.h> 6552f7358SJed Brown #include <exodusII.h> 7552f7358SJed Brown #endif 8552f7358SJed Brown 91e50132fSMatthew G. Knepley #include <petsc/private/viewerimpl.h> 101e50132fSMatthew G. Knepley #include <petsc/private/viewerexodusiiimpl.h> 11e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 12a1cb98faSBarry Smith /*@C 13a1cb98faSBarry Smith PETSC_VIEWER_EXODUSII_ - Creates an `PETSCVIEWEREXODUSII` `PetscViewer` shared by all processors in a communicator. 141e50132fSMatthew G. Knepley 1520f4b53cSBarry Smith Collective; No Fortran Support 161e50132fSMatthew G. Knepley 171e50132fSMatthew G. Knepley Input Parameter: 18a1cb98faSBarry Smith . comm - the MPI communicator to share the `PETSCVIEWEREXODUSII` `PetscViewer` 191e50132fSMatthew G. Knepley 201e50132fSMatthew G. Knepley Level: intermediate 211e50132fSMatthew G. Knepley 22a1cb98faSBarry Smith Note: 2320f4b53cSBarry Smith Unlike almost all other PETSc routines, `PETSC_VIEWER_EXODUSII_()` does not return 241e50132fSMatthew G. Knepley an error code. The GLVIS PetscViewer is usually used in the form 251e50132fSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm)); 261e50132fSMatthew G. Knepley 27a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewer`, `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()` 28a1cb98faSBarry Smith @*/ 29d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm) 30d71ae5a4SJacob Faibussowitsch { 311e50132fSMatthew G. Knepley PetscViewer viewer; 321e50132fSMatthew G. Knepley PetscErrorCode ierr; 331e50132fSMatthew G. Knepley 341e50132fSMatthew G. Knepley PetscFunctionBegin; 351e50132fSMatthew G. Knepley ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer); 369371c9d4SSatish Balay if (ierr) { 373ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 3811cc89d2SBarry Smith PetscFunctionReturn(NULL); 399371c9d4SSatish Balay } 401e50132fSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 419371c9d4SSatish Balay if (ierr) { 423ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 4311cc89d2SBarry Smith PetscFunctionReturn(NULL); 449371c9d4SSatish Balay } 451e50132fSMatthew G. Knepley PetscFunctionReturn(viewer); 461e50132fSMatthew G. Knepley } 471e50132fSMatthew G. Knepley 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer) 49d71ae5a4SJacob Faibussowitsch { 506823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data; 511e50132fSMatthew G. Knepley 521e50132fSMatthew G. Knepley PetscFunctionBegin; 539566063dSJacob Faibussowitsch if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename)); 549566063dSJacob Faibussowitsch if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid: %d\n", exo->exoid)); 559566063dSJacob Faibussowitsch if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode: %d\n", exo->btype)); 56197f6770SJed Brown if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order: %" PetscInt_FMT "\n", exo->order)); 573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 581e50132fSMatthew G. Knepley } 591e50132fSMatthew G. Knepley 60d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscViewer v, PetscOptionItems *PetscOptionsObject) 61d71ae5a4SJacob Faibussowitsch { 621e50132fSMatthew G. Knepley PetscFunctionBegin; 63d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options"); 64d0609cedSBarry Smith PetscOptionsHeadEnd(); 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 661e50132fSMatthew G. Knepley } 671e50132fSMatthew G. Knepley 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer) 69d71ae5a4SJacob Faibussowitsch { 701e50132fSMatthew G. Knepley PetscFunctionBegin; 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 721e50132fSMatthew G. Knepley } 731e50132fSMatthew G. Knepley 74d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer) 75d71ae5a4SJacob Faibussowitsch { 761e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 771e50132fSMatthew G. Knepley 781e50132fSMatthew G. Knepley PetscFunctionBegin; 7948a46eb9SPierre Jolivet if (exo->exoid >= 0) PetscCallExternal(ex_close, exo->exoid); 809566063dSJacob Faibussowitsch PetscCall(PetscFree(exo->filename)); 819566063dSJacob Faibussowitsch PetscCall(PetscFree(exo)); 829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 852e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetId_C", NULL)); 879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetOrder_C", NULL)); 889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerSetOrder_C", NULL)); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 901e50132fSMatthew G. Knepley } 911e50132fSMatthew G. Knepley 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[]) 93d71ae5a4SJacob Faibussowitsch { 941e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 951e50132fSMatthew G. Knepley PetscMPIInt rank; 961e50132fSMatthew G. Knepley int CPU_word_size, IO_word_size, EXO_mode; 976823f3c5SBlaise Bourdin MPI_Info mpi_info = MPI_INFO_NULL; 986823f3c5SBlaise Bourdin float EXO_version; 991e50132fSMatthew G. Knepley 1009566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 1011e50132fSMatthew G. Knepley CPU_word_size = sizeof(PetscReal); 1021e50132fSMatthew G. Knepley IO_word_size = sizeof(PetscReal); 1031e50132fSMatthew G. Knepley 1041e50132fSMatthew G. Knepley PetscFunctionBegin; 1056823f3c5SBlaise Bourdin if (exo->exoid >= 0) { 106792fecdfSBarry Smith PetscCallExternal(ex_close, exo->exoid); 1076823f3c5SBlaise Bourdin exo->exoid = -1; 1086823f3c5SBlaise Bourdin } 1099566063dSJacob Faibussowitsch if (exo->filename) PetscCall(PetscFree(exo->filename)); 1109566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &exo->filename)); 1111e50132fSMatthew G. Knepley switch (exo->btype) { 112d71ae5a4SJacob Faibussowitsch case FILE_MODE_READ: 113d71ae5a4SJacob Faibussowitsch EXO_mode = EX_READ; 114d71ae5a4SJacob Faibussowitsch break; 1151e50132fSMatthew G. Knepley case FILE_MODE_APPEND: 1166823f3c5SBlaise Bourdin case FILE_MODE_UPDATE: 1176823f3c5SBlaise Bourdin case FILE_MODE_APPEND_UPDATE: 1186823f3c5SBlaise Bourdin /* Will fail if the file does not already exist */ 1196823f3c5SBlaise Bourdin EXO_mode = EX_WRITE; 1201e50132fSMatthew G. Knepley break; 1211e50132fSMatthew G. Knepley case FILE_MODE_WRITE: 1226823f3c5SBlaise Bourdin /* 1236823f3c5SBlaise Bourdin exodus only allows writing geometry upon file creation, so we will let DMView create the file. 1246823f3c5SBlaise Bourdin */ 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1261e50132fSMatthew G. Knepley break; 127d71ae5a4SJacob Faibussowitsch default: 128d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 1291e50132fSMatthew G. Knepley } 1301e50132fSMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 1311e50132fSMatthew G. Knepley EXO_mode += EX_ALL_INT64_API; 1321e50132fSMatthew G. Knepley #endif 1336823f3c5SBlaise Bourdin exo->exoid = ex_open_par(name, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PETSC_COMM_WORLD, mpi_info); 13408401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name); 1353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1361e50132fSMatthew G. Knepley } 1371e50132fSMatthew G. Knepley 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name) 139d71ae5a4SJacob Faibussowitsch { 1401e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1411e50132fSMatthew G. Knepley 1421e50132fSMatthew G. Knepley PetscFunctionBegin; 1431e50132fSMatthew G. Knepley *name = exo->filename; 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1451e50132fSMatthew G. Knepley } 1461e50132fSMatthew G. Knepley 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type) 148d71ae5a4SJacob Faibussowitsch { 1491e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1501e50132fSMatthew G. Knepley 1511e50132fSMatthew G. Knepley PetscFunctionBegin; 1521e50132fSMatthew G. Knepley exo->btype = type; 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1541e50132fSMatthew G. Knepley } 1551e50132fSMatthew G. Knepley 156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type) 157d71ae5a4SJacob Faibussowitsch { 1581e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1591e50132fSMatthew G. Knepley 1601e50132fSMatthew G. Knepley PetscFunctionBegin; 1611e50132fSMatthew G. Knepley *type = exo->btype; 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1631e50132fSMatthew G. Knepley } 1641e50132fSMatthew G. Knepley 165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid) 166d71ae5a4SJacob Faibussowitsch { 1671e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1681e50132fSMatthew G. Knepley 1691e50132fSMatthew G. Knepley PetscFunctionBegin; 1701e50132fSMatthew G. Knepley *exoid = exo->exoid; 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1721e50132fSMatthew G. Knepley } 1731e50132fSMatthew G. Knepley 174d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order) 175d71ae5a4SJacob Faibussowitsch { 1766823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1771e50132fSMatthew G. Knepley 1781e50132fSMatthew G. Knepley PetscFunctionBegin; 1796823f3c5SBlaise Bourdin *order = exo->order; 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1816823f3c5SBlaise Bourdin } 1826823f3c5SBlaise Bourdin 183d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order) 184d71ae5a4SJacob Faibussowitsch { 1856823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 1866823f3c5SBlaise Bourdin 1876823f3c5SBlaise Bourdin PetscFunctionBegin; 1886823f3c5SBlaise Bourdin exo->order = order; 1893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1901e50132fSMatthew G. Knepley } 1911e50132fSMatthew G. Knepley 1921e50132fSMatthew G. Knepley /*MC 19300373969SVaclav Hapla PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file 1941e50132fSMatthew G. Knepley 19520f4b53cSBarry Smith Level: beginner 19620f4b53cSBarry Smith 197db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`, 198db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 1991e50132fSMatthew G. Knepley M*/ 2001e50132fSMatthew G. Knepley 201d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v) 202d71ae5a4SJacob Faibussowitsch { 2031e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo; 2041e50132fSMatthew G. Knepley 2051e50132fSMatthew G. Knepley PetscFunctionBegin; 2064dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&exo)); 2071e50132fSMatthew G. Knepley 2081e50132fSMatthew G. Knepley v->data = (void *)exo; 2091e50132fSMatthew G. Knepley v->ops->destroy = PetscViewerDestroy_ExodusII; 2101e50132fSMatthew G. Knepley v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII; 2111e50132fSMatthew G. Knepley v->ops->setup = PetscViewerSetUp_ExodusII; 2121e50132fSMatthew G. Knepley v->ops->view = PetscViewerView_ExodusII; 2131e50132fSMatthew G. Knepley v->ops->flush = 0; 2141e50132fSMatthew G. Knepley exo->btype = (PetscFileMode)-1; 2151e50132fSMatthew G. Knepley exo->filename = 0; 2161e50132fSMatthew G. Knepley exo->exoid = -1; 2171e50132fSMatthew G. Knepley 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_ExodusII)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_ExodusII)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_ExodusII)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_ExodusII)); 2229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetId_C", PetscViewerExodusIIGetId_ExodusII)); 2239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerSetOrder_C", PetscViewerExodusIISetOrder_ExodusII)); 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetOrder_C", PetscViewerExodusIIGetOrder_ExodusII)); 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2261e50132fSMatthew G. Knepley } 2271e50132fSMatthew G. Knepley 2281e50132fSMatthew G. Knepley /* 22978569bb4SMatthew G. Knepley EXOGetVarIndex - Locate a result in an exodus file based on its name 23078569bb4SMatthew G. Knepley 2316823f3c5SBlaise Bourdin Collective 23278569bb4SMatthew G. Knepley 23378569bb4SMatthew G. Knepley Input Parameters: 23478569bb4SMatthew G. Knepley + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 23578569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 23678569bb4SMatthew G. Knepley - name - the name of the result 23778569bb4SMatthew G. Knepley 2382fe279fdSBarry Smith Output Parameter: 23978569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result 24078569bb4SMatthew G. Knepley 24120f4b53cSBarry Smith Level: beginner 24220f4b53cSBarry Smith 24378569bb4SMatthew G. Knepley Notes: 24478569bb4SMatthew G. Knepley The exodus variable index is obtained by comparing name and the 24578569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance if name is "V" 24678569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 24778569bb4SMatthew G. Knepley amongst all variables of type obj_type. 24878569bb4SMatthew G. Knepley 249c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 25078569bb4SMatthew G. Knepley */ 251d71ae5a4SJacob Faibussowitsch PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 252d71ae5a4SJacob Faibussowitsch { 2536de834b4SMatthew G. Knepley int num_vars, i, j; 25478569bb4SMatthew G. Knepley char ext_name[MAX_STR_LENGTH + 1], var_name[MAX_STR_LENGTH + 1]; 25578569bb4SMatthew G. Knepley const int num_suffix = 5; 25678569bb4SMatthew G. Knepley char *suffix[5]; 25778569bb4SMatthew G. Knepley PetscBool flg; 25878569bb4SMatthew G. Knepley 25978569bb4SMatthew G. Knepley PetscFunctionBegin; 26078569bb4SMatthew G. Knepley suffix[0] = (char *)""; 26178569bb4SMatthew G. Knepley suffix[1] = (char *)"_X"; 26278569bb4SMatthew G. Knepley suffix[2] = (char *)"_XX"; 26378569bb4SMatthew G. Knepley suffix[3] = (char *)"_1"; 26478569bb4SMatthew G. Knepley suffix[4] = (char *)"_11"; 26578569bb4SMatthew G. Knepley 2666823f3c5SBlaise Bourdin *varIndex = -1; 267792fecdfSBarry Smith PetscCallExternal(ex_get_variable_param, exoid, obj_type, &num_vars); 26878569bb4SMatthew G. Knepley for (i = 0; i < num_vars; ++i) { 269792fecdfSBarry Smith PetscCallExternal(ex_get_variable_name, exoid, obj_type, i + 1, var_name); 27078569bb4SMatthew G. Knepley for (j = 0; j < num_suffix; ++j) { 2719566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH)); 2729566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH)); 2739566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(ext_name, var_name, &flg)); 274ad540459SPierre Jolivet if (flg) *varIndex = i + 1; 2756823f3c5SBlaise Bourdin } 2766823f3c5SBlaise Bourdin } 2773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27878569bb4SMatthew G. Knepley } 27978569bb4SMatthew G. Knepley 28078569bb4SMatthew G. Knepley /* 28120f4b53cSBarry Smith DMView_PlexExodusII - Write a `DM` to disk in exodus format 28278569bb4SMatthew G. Knepley 28320f4b53cSBarry Smith Collective 28478569bb4SMatthew G. Knepley 28578569bb4SMatthew G. Knepley Input Parameters: 28678569bb4SMatthew G. Knepley + dm - The dm to be written 2876823f3c5SBlaise Bourdin . viewer - an exodusII viewer 28878569bb4SMatthew G. Knepley 28920f4b53cSBarry Smith Level: beginner 29020f4b53cSBarry Smith 29178569bb4SMatthew G. Knepley Notes: 29278569bb4SMatthew G. Knepley Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 29378569bb4SMatthew G. Knepley consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 29478569bb4SMatthew G. Knepley 29520f4b53cSBarry Smith If `dm` has been distributed, only the part of the `DM` on MPI rank 0 (including "ghost" cells and vertices) 2966823f3c5SBlaise Bourdin will be written. 29778569bb4SMatthew G. Knepley 29820f4b53cSBarry Smith `DMPLEX` only represents geometry while most post-processing software expect that a mesh also provides information 29978569bb4SMatthew G. Knepley on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 30020f4b53cSBarry Smith The order of the mesh shall be set using `PetscViewerExodusIISetOrder()` 30120f4b53cSBarry Smith It should be extended to use `PetscFE` objects. 30278569bb4SMatthew G. Knepley 3036823f3c5SBlaise Bourdin This function will only handle TRI, TET, QUAD, and HEX cells. 30478569bb4SMatthew G. Knepley 30520f4b53cSBarry Smith .seealso: `DMPLEX` 30678569bb4SMatthew G. Knepley */ 307d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer) 308d71ae5a4SJacob Faibussowitsch { 3099371c9d4SSatish Balay enum ElemType { 3109371c9d4SSatish Balay TRI, 3119371c9d4SSatish Balay QUAD, 3129371c9d4SSatish Balay TET, 3139371c9d4SSatish Balay HEX 3149371c9d4SSatish Balay }; 31578569bb4SMatthew G. Knepley MPI_Comm comm; 3166823f3c5SBlaise Bourdin PetscInt degree; /* the order of the mesh */ 31778569bb4SMatthew G. Knepley /* Connectivity Variables */ 31878569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 31978569bb4SMatthew G. Knepley /* Cell Sets */ 32078569bb4SMatthew G. Knepley DMLabel csLabel; 32178569bb4SMatthew G. Knepley IS csIS; 32278569bb4SMatthew G. Knepley const PetscInt *csIdx; 32378569bb4SMatthew G. Knepley PetscInt num_cs, cs; 32478569bb4SMatthew G. Knepley enum ElemType *type; 325fe209ef9SBlaise Bourdin PetscBool hasLabel; 32678569bb4SMatthew G. Knepley /* Coordinate Variables */ 32778569bb4SMatthew G. Knepley DM cdm; 3286823f3c5SBlaise Bourdin PetscSection coordSection; 32978569bb4SMatthew G. Knepley Vec coord; 33078569bb4SMatthew G. Knepley PetscInt **nodes; 331eae8a387SMatthew G. Knepley PetscInt depth, d, dim, skipCells = 0; 33278569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 33378569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 33478569bb4SMatthew G. Knepley PetscMPIInt rank, size; 33578569bb4SMatthew G. Knepley const char *dmName; 336fe209ef9SBlaise Bourdin PetscInt nodesTriP1[4] = {3, 0, 0, 0}; 337fe209ef9SBlaise Bourdin PetscInt nodesTriP2[4] = {3, 3, 0, 0}; 338fe209ef9SBlaise Bourdin PetscInt nodesQuadP1[4] = {4, 0, 0, 0}; 339fe209ef9SBlaise Bourdin PetscInt nodesQuadP2[4] = {4, 4, 0, 1}; 340fe209ef9SBlaise Bourdin PetscInt nodesTetP1[4] = {4, 0, 0, 0}; 341fe209ef9SBlaise Bourdin PetscInt nodesTetP2[4] = {4, 6, 0, 0}; 342fe209ef9SBlaise Bourdin PetscInt nodesHexP1[4] = {8, 0, 0, 0}; 343fe209ef9SBlaise Bourdin PetscInt nodesHexP2[4] = {8, 12, 6, 1}; 3446823f3c5SBlaise Bourdin int CPU_word_size, IO_word_size, EXO_mode; 3456823f3c5SBlaise Bourdin float EXO_version; 3466823f3c5SBlaise Bourdin 3476823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 3486823f3c5SBlaise Bourdin 34978569bb4SMatthew G. Knepley PetscFunctionBegin; 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3536823f3c5SBlaise Bourdin 3546823f3c5SBlaise Bourdin /* 3556823f3c5SBlaise Bourdin Creating coordSection is a collective operation so we do it somewhat out of sequence 3566823f3c5SBlaise Bourdin */ 3579566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &coordSection)); 3589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(dm)); 359bdbfaa5dSBlaise Bourdin /* 360bdbfaa5dSBlaise Bourdin Check that all points are on rank 0 since we don't know how to save distributed DM in exodus format 361bdbfaa5dSBlaise Bourdin */ 362bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 363bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 364bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 365bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 366bdbfaa5dSBlaise Bourdin numCells = cEnd - cStart; 367bdbfaa5dSBlaise Bourdin numEdges = eEnd - eStart; 368bdbfaa5dSBlaise Bourdin numVertices = vEnd - vStart; 369bdbfaa5dSBlaise Bourdin PetscCheck(!(rank && (numCells || numEdges || numVertices)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Writing distributed DM in exodusII format not supported"); 370c5853193SPierre Jolivet if (rank == 0) { 3716823f3c5SBlaise Bourdin switch (exo->btype) { 3726823f3c5SBlaise Bourdin case FILE_MODE_READ: 3736823f3c5SBlaise Bourdin case FILE_MODE_APPEND: 3746823f3c5SBlaise Bourdin case FILE_MODE_UPDATE: 3756823f3c5SBlaise Bourdin case FILE_MODE_APPEND_UPDATE: 3766823f3c5SBlaise Bourdin /* exodusII does not allow writing geometry to an existing file */ 37798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename); 3786823f3c5SBlaise Bourdin case FILE_MODE_WRITE: 3796823f3c5SBlaise Bourdin /* Create an empty file if one already exists*/ 3806823f3c5SBlaise Bourdin EXO_mode = EX_CLOBBER; 3816823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 3826823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 3836823f3c5SBlaise Bourdin #endif 3846823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 3856823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 3866823f3c5SBlaise Bourdin exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size); 38708401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename); 3886823f3c5SBlaise Bourdin 3896823f3c5SBlaise Bourdin break; 390d71ae5a4SJacob Faibussowitsch default: 391d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3926823f3c5SBlaise Bourdin } 3936823f3c5SBlaise Bourdin 39478569bb4SMatthew G. Knepley /* --- Get DM info --- */ 3959566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &dmName)); 3969566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 3979566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3989566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 3999371c9d4SSatish Balay if (depth == 3) { 4009371c9d4SSatish Balay numFaces = fEnd - fStart; 4019371c9d4SSatish Balay } else { 4029371c9d4SSatish Balay numFaces = 0; 4039371c9d4SSatish Balay } 4049566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs)); 4059566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs)); 4069566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs)); 4079566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coord)); 4089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 40978569bb4SMatthew G. Knepley if (num_cs > 0) { 4109566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel)); 4119566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(csLabel, &csIS)); 4129566063dSJacob Faibussowitsch PetscCall(ISGetIndices(csIS, &csIdx)); 41378569bb4SMatthew G. Knepley } 4149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &nodes)); 41578569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 4169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &type)); 41778569bb4SMatthew G. Knepley numNodes = numVertices; 4186823f3c5SBlaise Bourdin 4199566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetOrder(viewer, °ree)); 420ad540459SPierre Jolivet if (degree == 2) numNodes += numEdges; 42178569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 42278569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 42378569bb4SMatthew G. Knepley IS stratumIS; 42478569bb4SMatthew G. Knepley const PetscInt *cells; 42578569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 42678569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 42778569bb4SMatthew G. Knepley 4289566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 4299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 4309566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 4319566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 43278569bb4SMatthew G. Knepley switch (dim) { 43378569bb4SMatthew G. Knepley case 2: 4349371c9d4SSatish Balay if (closureSize == 3 * dim) { 4359371c9d4SSatish Balay type[cs] = TRI; 4369371c9d4SSatish Balay } else if (closureSize == 4 * dim) { 4379371c9d4SSatish Balay type[cs] = QUAD; 4389371c9d4SSatish Balay } 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); 43978569bb4SMatthew G. Knepley break; 44078569bb4SMatthew G. Knepley case 3: 4419371c9d4SSatish Balay if (closureSize == 4 * dim) { 4429371c9d4SSatish Balay type[cs] = TET; 4439371c9d4SSatish Balay } else if (closureSize == 8 * dim) { 4449371c9d4SSatish Balay type[cs] = HEX; 4459371c9d4SSatish Balay } 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); 44678569bb4SMatthew G. Knepley break; 447d71ae5a4SJacob Faibussowitsch default: 448d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim); 44978569bb4SMatthew G. Knepley } 450ad540459SPierre Jolivet if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize; 4519371c9d4SSatish Balay if ((degree == 2) && (type[cs] == HEX)) { 4529371c9d4SSatish Balay numNodes += csSize; 4539371c9d4SSatish Balay numNodes += numFaces; 4549371c9d4SSatish Balay } 4559566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 45678569bb4SMatthew G. Knepley /* Set nodes and Element type */ 45778569bb4SMatthew G. Knepley if (type[cs] == TRI) { 45878569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 45978569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 46078569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 46178569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 46278569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 46378569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 46478569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 46578569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 46678569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 46778569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 46878569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 46978569bb4SMatthew G. Knepley } 47078569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 47178569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3] * csSize; 47278569bb4SMatthew G. Knepley 4739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 4749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 47578569bb4SMatthew G. Knepley } 476bdbfaa5dSBlaise Bourdin if (num_cs) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs); 47778569bb4SMatthew G. Knepley /* --- Connectivity --- */ 47878569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 47978569bb4SMatthew G. Knepley IS stratumIS; 48078569bb4SMatthew G. Knepley const PetscInt *cells; 481fe209ef9SBlaise Bourdin PetscInt *connect, off = 0; 482eae8a387SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 48378569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 48478569bb4SMatthew G. Knepley char *elem_type = NULL; 48578569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 48678569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 48778569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 48878569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 48978569bb4SMatthew G. Knepley 4909566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 4919566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 4929566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 49378569bb4SMatthew G. Knepley /* Set Element type */ 49478569bb4SMatthew G. Knepley if (type[cs] == TRI) { 49578569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 49678569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 49778569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 49878569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 49978569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 50078569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 50178569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 50278569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 50378569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 50478569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 50578569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 50678569bb4SMatthew G. Knepley } 50778569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 5089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect)); 509792fecdfSBarry Smith PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1); 51078569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 51178569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 51278569bb4SMatthew G. Knepley if (depth > 1) { 51378569bb4SMatthew G. Knepley if (dim == 2) { 5149566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure)); 51578569bb4SMatthew G. Knepley } else if (dim == 3) { 51678569bb4SMatthew G. Knepley PetscInt *closure = NULL; 51778569bb4SMatthew G. Knepley 5189566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure)); 5199566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 52078569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 5219566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 52278569bb4SMatthew G. Knepley } 52378569bb4SMatthew G. Knepley } 52478569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 52578569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 52678569bb4SMatthew G. Knepley PetscInt *closure = NULL; 52778569bb4SMatthew G. Knepley PetscInt temp, i; 52878569bb4SMatthew G. Knepley 5299566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 53078569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 53178569bb4SMatthew G. Knepley if (i < nodes[cs][0]) { /* Vertices */ 532fe209ef9SBlaise Bourdin connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1; 533fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 53478569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */ 535fe209ef9SBlaise Bourdin connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1; 536fe209ef9SBlaise Bourdin if (nodes[cs][2] == 0) connect[i + off] -= numFaces; 537fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 53878569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */ 539fe209ef9SBlaise Bourdin connect[i + off] = closure[0] + 1; 540fe209ef9SBlaise Bourdin connect[i + off] -= skipCells; 54178569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */ 542fe209ef9SBlaise Bourdin connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1; 543fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 54478569bb4SMatthew G. Knepley } else { 545fe209ef9SBlaise Bourdin connect[i + off] = -1; 54678569bb4SMatthew G. Knepley } 54778569bb4SMatthew G. Knepley } 54878569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 54978569bb4SMatthew G. Knepley if (type[cs] == TET) { 5509371c9d4SSatish Balay temp = connect[0 + off]; 5519371c9d4SSatish Balay connect[0 + off] = connect[1 + off]; 5529371c9d4SSatish Balay connect[1 + off] = temp; 55378569bb4SMatthew G. Knepley if (degree == 2) { 5549371c9d4SSatish Balay temp = connect[5 + off]; 5559371c9d4SSatish Balay connect[5 + off] = connect[6 + off]; 5569371c9d4SSatish Balay connect[6 + off] = temp; 5579371c9d4SSatish Balay temp = connect[7 + off]; 5589371c9d4SSatish Balay connect[7 + off] = connect[8 + off]; 5599371c9d4SSatish Balay connect[8 + off] = temp; 56078569bb4SMatthew G. Knepley } 56178569bb4SMatthew G. Knepley } 56278569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 56378569bb4SMatthew G. Knepley if (type[cs] == HEX) { 5649371c9d4SSatish Balay temp = connect[1 + off]; 5659371c9d4SSatish Balay connect[1 + off] = connect[3 + off]; 5669371c9d4SSatish Balay connect[3 + off] = temp; 56778569bb4SMatthew G. Knepley if (degree == 2) { 5689371c9d4SSatish Balay temp = connect[8 + off]; 5699371c9d4SSatish Balay connect[8 + off] = connect[11 + off]; 5709371c9d4SSatish Balay connect[11 + off] = temp; 5719371c9d4SSatish Balay temp = connect[9 + off]; 5729371c9d4SSatish Balay connect[9 + off] = connect[10 + off]; 5739371c9d4SSatish Balay connect[10 + off] = temp; 5749371c9d4SSatish Balay temp = connect[16 + off]; 5759371c9d4SSatish Balay connect[16 + off] = connect[17 + off]; 5769371c9d4SSatish Balay connect[17 + off] = temp; 5779371c9d4SSatish Balay temp = connect[18 + off]; 5789371c9d4SSatish Balay connect[18 + off] = connect[19 + off]; 5799371c9d4SSatish Balay connect[19 + off] = temp; 58078569bb4SMatthew G. Knepley 5819371c9d4SSatish Balay temp = connect[12 + off]; 5829371c9d4SSatish Balay connect[12 + off] = connect[16 + off]; 5839371c9d4SSatish Balay connect[16 + off] = temp; 5849371c9d4SSatish Balay temp = connect[13 + off]; 5859371c9d4SSatish Balay connect[13 + off] = connect[17 + off]; 5869371c9d4SSatish Balay connect[17 + off] = temp; 5879371c9d4SSatish Balay temp = connect[14 + off]; 5889371c9d4SSatish Balay connect[14 + off] = connect[18 + off]; 5899371c9d4SSatish Balay connect[18 + off] = temp; 5909371c9d4SSatish Balay temp = connect[15 + off]; 5919371c9d4SSatish Balay connect[15 + off] = connect[19 + off]; 5929371c9d4SSatish Balay connect[19 + off] = temp; 59378569bb4SMatthew G. Knepley 5949371c9d4SSatish Balay temp = connect[23 + off]; 5959371c9d4SSatish Balay connect[23 + off] = connect[26 + off]; 5969371c9d4SSatish Balay connect[26 + off] = temp; 5979371c9d4SSatish Balay temp = connect[24 + off]; 5989371c9d4SSatish Balay connect[24 + off] = connect[25 + off]; 5999371c9d4SSatish Balay connect[25 + off] = temp; 6009371c9d4SSatish Balay temp = connect[25 + off]; 6019371c9d4SSatish Balay connect[25 + off] = connect[26 + off]; 6029371c9d4SSatish Balay connect[26 + off] = temp; 60378569bb4SMatthew G. Knepley } 60478569bb4SMatthew G. Knepley } 605fe209ef9SBlaise Bourdin off += connectSize; 6069566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 60778569bb4SMatthew G. Knepley } 608792fecdfSBarry Smith PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0); 60978569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0) * csSize; 6109566063dSJacob Faibussowitsch PetscCall(PetscFree(connect)); 6119566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 6129566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 61378569bb4SMatthew G. Knepley } 6149566063dSJacob Faibussowitsch PetscCall(PetscFree(type)); 61578569bb4SMatthew G. Knepley /* --- Coordinates --- */ 6169566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd)); 6171e50132fSMatthew G. Knepley if (num_cs) { 61878569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 6199566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 62048a46eb9SPierre Jolivet for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0)); 62178569bb4SMatthew G. Knepley } 6221e50132fSMatthew G. Knepley } 62378569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 62478569bb4SMatthew G. Knepley IS stratumIS; 62578569bb4SMatthew G. Knepley const PetscInt *cells; 62678569bb4SMatthew G. Knepley PetscInt csSize, c; 62778569bb4SMatthew G. Knepley 6289566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 6299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 6309566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 63148a46eb9SPierre Jolivet for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0)); 6329566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 6339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 63478569bb4SMatthew G. Knepley } 635bdbfaa5dSBlaise Bourdin if (num_cs) { 6369566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(csIS, &csIdx)); 6379566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS)); 63878569bb4SMatthew G. Knepley } 6399566063dSJacob Faibussowitsch PetscCall(PetscFree(nodes)); 6409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 641bdbfaa5dSBlaise Bourdin if (numNodes) { 64278569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 643233c95e0SFrancesco Ballarin PetscScalar *closure, *cval; 644233c95e0SFrancesco Ballarin PetscReal *coords; 64578569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 64678569bb4SMatthew G. Knepley 6476823f3c5SBlaise Bourdin /* There can't be more than 24 values in the closure of a point for the coord coordSection */ 6489566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure)); 6499566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord)); 6509566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 65178569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 6529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, p, &hasDof)); 65378569bb4SMatthew G. Knepley if (hasDof) { 65478569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 65578569bb4SMatthew G. Knepley 6569566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure)); 65778569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 65878569bb4SMatthew G. Knepley cval[d] = 0.0; 65978569bb4SMatthew G. Knepley for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d]; 660233c95e0SFrancesco Ballarin coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize; 66178569bb4SMatthew G. Knepley } 66278569bb4SMatthew G. Knepley ++n; 66378569bb4SMatthew G. Knepley } 66478569bb4SMatthew G. Knepley } 665792fecdfSBarry Smith PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]); 6669566063dSJacob Faibussowitsch PetscCall(PetscFree3(coords, cval, closure)); 667792fecdfSBarry Smith PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames); 66878569bb4SMatthew G. Knepley } 6696823f3c5SBlaise Bourdin 670fe209ef9SBlaise Bourdin /* --- Node Sets/Vertex Sets --- */ 6713ba16761SJacob Faibussowitsch PetscCall(DMHasLabel(dm, "Vertex Sets", &hasLabel)); 672fe209ef9SBlaise Bourdin if (hasLabel) { 673fe209ef9SBlaise Bourdin PetscInt i, vs, vsSize; 674fe209ef9SBlaise Bourdin const PetscInt *vsIdx, *vertices; 675fe209ef9SBlaise Bourdin PetscInt *nodeList; 676fe209ef9SBlaise Bourdin IS vsIS, stratumIS; 677fe209ef9SBlaise Bourdin DMLabel vsLabel; 6789566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel)); 6799566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(vsLabel, &vsIS)); 6809566063dSJacob Faibussowitsch PetscCall(ISGetIndices(vsIS, &vsIdx)); 681fe209ef9SBlaise Bourdin for (vs = 0; vs < num_vs; ++vs) { 6829566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS)); 6839566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &vertices)); 6849566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &vsSize)); 6859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(vsSize, &nodeList)); 686ad540459SPierre Jolivet for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1; 687792fecdfSBarry Smith PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0); 688792fecdfSBarry Smith PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL); 6899566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &vertices)); 6909566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 6919566063dSJacob Faibussowitsch PetscCall(PetscFree(nodeList)); 692fe209ef9SBlaise Bourdin } 6939566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(vsIS, &vsIdx)); 6949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vsIS)); 695fe209ef9SBlaise Bourdin } 696fe209ef9SBlaise Bourdin /* --- Side Sets/Face Sets --- */ 6979566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel)); 698fe209ef9SBlaise Bourdin if (hasLabel) { 699fe209ef9SBlaise Bourdin PetscInt i, j, fs, fsSize; 700fe209ef9SBlaise Bourdin const PetscInt *fsIdx, *faces; 701fe209ef9SBlaise Bourdin IS fsIS, stratumIS; 702fe209ef9SBlaise Bourdin DMLabel fsLabel; 703fe209ef9SBlaise Bourdin PetscInt numPoints, *points; 704fe209ef9SBlaise Bourdin PetscInt elem_list_size = 0; 705fe209ef9SBlaise Bourdin PetscInt *elem_list, *elem_ind, *side_list; 706fe209ef9SBlaise Bourdin 7079566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel)); 708fe209ef9SBlaise Bourdin /* Compute size of Node List and Element List */ 7099566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(fsLabel, &fsIS)); 7109566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fsIS, &fsIdx)); 711fe209ef9SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) { 7129566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 7139566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 714fe209ef9SBlaise Bourdin elem_list_size += fsSize; 7159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 716fe209ef9SBlaise Bourdin } 717bdbfaa5dSBlaise Bourdin if (num_fs) { 718d0609cedSBarry Smith PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list)); 719fe209ef9SBlaise Bourdin elem_ind[0] = 0; 720fe209ef9SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) { 7219566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 7229566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &faces)); 7239566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 724fe209ef9SBlaise Bourdin /* Set Parameters */ 725792fecdfSBarry Smith PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0); 726fe209ef9SBlaise Bourdin /* Indices */ 727ad540459SPierre Jolivet if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize; 728fe209ef9SBlaise Bourdin 729fe209ef9SBlaise Bourdin for (i = 0; i < fsSize; ++i) { 730fe209ef9SBlaise Bourdin /* Element List */ 731fe209ef9SBlaise Bourdin points = NULL; 7329566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 733fe209ef9SBlaise Bourdin elem_list[elem_ind[fs] + i] = points[2] + 1; 7349566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 735fe209ef9SBlaise Bourdin 736fe209ef9SBlaise Bourdin /* Side List */ 737fe209ef9SBlaise Bourdin points = NULL; 7389566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 739fe209ef9SBlaise Bourdin for (j = 1; j < numPoints; ++j) { 740ad540459SPierre Jolivet if (points[j * 2] == faces[i]) break; 741fe209ef9SBlaise Bourdin } 742fe209ef9SBlaise Bourdin /* Convert HEX sides */ 743fe209ef9SBlaise Bourdin if (numPoints == 27) { 7449371c9d4SSatish Balay if (j == 1) { 7459371c9d4SSatish Balay j = 5; 7469371c9d4SSatish Balay } else if (j == 2) { 7479371c9d4SSatish Balay j = 6; 7489371c9d4SSatish Balay } else if (j == 3) { 7499371c9d4SSatish Balay j = 1; 7509371c9d4SSatish Balay } else if (j == 4) { 7519371c9d4SSatish Balay j = 3; 7529371c9d4SSatish Balay } else if (j == 5) { 7539371c9d4SSatish Balay j = 2; 7549371c9d4SSatish Balay } else if (j == 6) { 7559371c9d4SSatish Balay j = 4; 7569371c9d4SSatish Balay } 757fe209ef9SBlaise Bourdin } 758fe209ef9SBlaise Bourdin /* Convert TET sides */ 759fe209ef9SBlaise Bourdin if (numPoints == 15) { 760fe209ef9SBlaise Bourdin --j; 761ad540459SPierre Jolivet if (j == 0) j = 4; 762fe209ef9SBlaise Bourdin } 763fe209ef9SBlaise Bourdin side_list[elem_ind[fs] + i] = j; 7649566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 765fe209ef9SBlaise Bourdin } 7669566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &faces)); 7679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 768fe209ef9SBlaise Bourdin } 7699566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fsIS, &fsIdx)); 7709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fsIS)); 771fe209ef9SBlaise Bourdin 772fe209ef9SBlaise Bourdin /* Put side sets */ 77348a46eb9SPierre Jolivet 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]]); 7749566063dSJacob Faibussowitsch PetscCall(PetscFree3(elem_ind, elem_list, side_list)); 775fe209ef9SBlaise Bourdin } 776bdbfaa5dSBlaise Bourdin } 7776823f3c5SBlaise Bourdin /* 7786823f3c5SBlaise Bourdin close the exodus file 7796823f3c5SBlaise Bourdin */ 7806823f3c5SBlaise Bourdin ex_close(exo->exoid); 7816823f3c5SBlaise Bourdin exo->exoid = -1; 7826823f3c5SBlaise Bourdin } 7839566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&coordSection)); 7846823f3c5SBlaise Bourdin 7856823f3c5SBlaise Bourdin /* 7866823f3c5SBlaise Bourdin reopen the file in parallel 7876823f3c5SBlaise Bourdin */ 7886823f3c5SBlaise Bourdin EXO_mode = EX_WRITE; 7896823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 7906823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 7916823f3c5SBlaise Bourdin #endif 7926823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 7936823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 7946823f3c5SBlaise Bourdin exo->exoid = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL); 79508401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename); 7963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7976823f3c5SBlaise Bourdin } 7986823f3c5SBlaise Bourdin 7996823f3c5SBlaise Bourdin /* 8006823f3c5SBlaise Bourdin VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 8016823f3c5SBlaise Bourdin 80220f4b53cSBarry Smith Collective 8036823f3c5SBlaise Bourdin 8046823f3c5SBlaise Bourdin Input Parameters: 8056823f3c5SBlaise Bourdin + v - The vector to be written 8066823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 8076823f3c5SBlaise Bourdin 80820f4b53cSBarry Smith Level: beginner 80920f4b53cSBarry Smith 8106823f3c5SBlaise Bourdin Notes: 8116823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 8126823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 8136823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 8146823f3c5SBlaise Bourdin amongst all variables. 8156823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 8166823f3c5SBlaise Bourdin possibly corrupting the file 8176823f3c5SBlaise Bourdin 818c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 8196823f3c5SBlaise Bourdin @*/ 820d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer) 821d71ae5a4SJacob Faibussowitsch { 8226823f3c5SBlaise Bourdin DM dm; 8236823f3c5SBlaise Bourdin MPI_Comm comm; 8246823f3c5SBlaise Bourdin PetscMPIInt rank; 8256823f3c5SBlaise Bourdin 8266823f3c5SBlaise Bourdin int exoid, offsetN = 0, offsetZ = 0; 8276823f3c5SBlaise Bourdin const char *vecname; 8286823f3c5SBlaise Bourdin PetscInt step; 8296823f3c5SBlaise Bourdin 8306823f3c5SBlaise Bourdin PetscFunctionBegin; 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8339566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 8349566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 8366823f3c5SBlaise Bourdin 8379566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 8389566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 8399566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 8401dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 8416823f3c5SBlaise Bourdin if (offsetN > 0) { 8429566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 8436823f3c5SBlaise Bourdin } else if (offsetZ > 0) { 8449566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 84598921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8476823f3c5SBlaise Bourdin } 8486823f3c5SBlaise Bourdin 8496823f3c5SBlaise Bourdin /* 8506823f3c5SBlaise Bourdin VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 8516823f3c5SBlaise Bourdin 85220f4b53cSBarry Smith Collective 8536823f3c5SBlaise Bourdin 8546823f3c5SBlaise Bourdin Input Parameters: 8556823f3c5SBlaise Bourdin + v - The vector to be written 8566823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 8576823f3c5SBlaise Bourdin 85820f4b53cSBarry Smith Level: beginner 85920f4b53cSBarry Smith 8606823f3c5SBlaise Bourdin Notes: 8616823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 8626823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 8636823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 8646823f3c5SBlaise Bourdin amongst all variables. 8656823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 8666823f3c5SBlaise Bourdin possibly corrupting the file 8676823f3c5SBlaise Bourdin 868c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 8696823f3c5SBlaise Bourdin @*/ 870d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer) 871d71ae5a4SJacob Faibussowitsch { 8726823f3c5SBlaise Bourdin DM dm; 8736823f3c5SBlaise Bourdin MPI_Comm comm; 8746823f3c5SBlaise Bourdin PetscMPIInt rank; 8756823f3c5SBlaise Bourdin 8766823f3c5SBlaise Bourdin int exoid, offsetN = 0, offsetZ = 0; 8776823f3c5SBlaise Bourdin const char *vecname; 8786823f3c5SBlaise Bourdin PetscInt step; 8796823f3c5SBlaise Bourdin 8806823f3c5SBlaise Bourdin PetscFunctionBegin; 8819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 8829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8839566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 8849566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 8866823f3c5SBlaise Bourdin 8879566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 8889566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 8899566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 8901dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 8911dca8a05SBarry Smith if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 8921dca8a05SBarry Smith else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 8931dca8a05SBarry Smith else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 8943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89578569bb4SMatthew G. Knepley } 89678569bb4SMatthew G. Knepley 89778569bb4SMatthew G. Knepley /* 89878569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 89978569bb4SMatthew G. Knepley 90020f4b53cSBarry Smith Collective 90178569bb4SMatthew G. Knepley 90278569bb4SMatthew G. Knepley Input Parameters: 90378569bb4SMatthew G. Knepley + v - The vector to be written 90478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9056823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 9066823f3c5SBlaise Bourdin - offset - the location of the variable in the file 90778569bb4SMatthew G. Knepley 90820f4b53cSBarry Smith Level: beginner 90920f4b53cSBarry Smith 91078569bb4SMatthew G. Knepley Notes: 91178569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 91278569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 91378569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 91478569bb4SMatthew G. Knepley amongst all nodal variables. 91578569bb4SMatthew G. Knepley 916c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 91778569bb4SMatthew G. Knepley @*/ 918d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 919d71ae5a4SJacob Faibussowitsch { 92078569bb4SMatthew G. Knepley MPI_Comm comm; 92178569bb4SMatthew G. Knepley PetscMPIInt size; 92278569bb4SMatthew G. Knepley DM dm; 92378569bb4SMatthew G. Knepley Vec vNatural, vComp; 92422a7544eSVaclav Hapla const PetscScalar *varray; 92578569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 92678569bb4SMatthew G. Knepley PetscBool useNatural; 92778569bb4SMatthew G. Knepley 92878569bb4SMatthew G. Knepley PetscFunctionBegin; 9299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 9309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 9319566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 9329566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 93378569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 93478569bb4SMatthew G. Knepley if (useNatural) { 935b7352c5cSAlexis Marboeuf PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 9369566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 9379566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 93878569bb4SMatthew G. Knepley } else { 93978569bb4SMatthew G. Knepley vNatural = v; 94078569bb4SMatthew G. Knepley } 9416823f3c5SBlaise Bourdin 94278569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 94378569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 94478569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 9459566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 9469566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 94778569bb4SMatthew G. Knepley if (bs == 1) { 9489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 949792fecdfSBarry Smith PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 9509566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 95178569bb4SMatthew G. Knepley } else { 95278569bb4SMatthew G. Knepley IS compIS; 95378569bb4SMatthew G. Knepley PetscInt c; 95478569bb4SMatthew G. Knepley 9559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 95678569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 9579566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 9589566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 9599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 960792fecdfSBarry Smith PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 9619566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 9629566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 96378569bb4SMatthew G. Knepley } 9649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 96578569bb4SMatthew G. Knepley } 966b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(VecDestroy(&vNatural)); 9673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96878569bb4SMatthew G. Knepley } 96978569bb4SMatthew G. Knepley 97078569bb4SMatthew G. Knepley /* 97178569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 97278569bb4SMatthew G. Knepley 97320f4b53cSBarry Smith Collective 97478569bb4SMatthew G. Knepley 97578569bb4SMatthew G. Knepley Input Parameters: 97678569bb4SMatthew G. Knepley + v - The vector to be written 97778569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9786823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 9796823f3c5SBlaise Bourdin - offset - the location of the variable in the file 98078569bb4SMatthew G. Knepley 98120f4b53cSBarry Smith Level: beginner 98220f4b53cSBarry Smith 98378569bb4SMatthew G. Knepley Notes: 98478569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 98578569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 98678569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 98778569bb4SMatthew G. Knepley amongst all nodal variables. 98878569bb4SMatthew G. Knepley 989db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 99078569bb4SMatthew G. Knepley */ 991d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 992d71ae5a4SJacob Faibussowitsch { 99378569bb4SMatthew G. Knepley MPI_Comm comm; 99478569bb4SMatthew G. Knepley PetscMPIInt size; 99578569bb4SMatthew G. Knepley DM dm; 99678569bb4SMatthew G. Knepley Vec vNatural, vComp; 99722a7544eSVaclav Hapla PetscScalar *varray; 99878569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 99978569bb4SMatthew G. Knepley PetscBool useNatural; 100078569bb4SMatthew G. Knepley 100178569bb4SMatthew G. Knepley PetscFunctionBegin; 10029566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 10039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10049566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 10059566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 100678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1007b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1008ad540459SPierre Jolivet else vNatural = v; 10096823f3c5SBlaise Bourdin 101078569bb4SMatthew G. Knepley /* Read local chunk from the file */ 10119566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 10129566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 101378569bb4SMatthew G. Knepley if (bs == 1) { 10149566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 1015792fecdfSBarry Smith PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 10169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 101778569bb4SMatthew G. Knepley } else { 101878569bb4SMatthew G. Knepley IS compIS; 101978569bb4SMatthew G. Knepley PetscInt c; 102078569bb4SMatthew G. Knepley 10219566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 102278569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 10239566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 10249566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 10259566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 1026792fecdfSBarry Smith PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 10279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 10289566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 102978569bb4SMatthew G. Knepley } 10309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 103178569bb4SMatthew G. Knepley } 103278569bb4SMatthew G. Knepley if (useNatural) { 10339566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 10349566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1035b7352c5cSAlexis Marboeuf PetscCall(VecDestroy(&vNatural)); 103678569bb4SMatthew G. Knepley } 10373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 103878569bb4SMatthew G. Knepley } 103978569bb4SMatthew G. Knepley 104078569bb4SMatthew G. Knepley /* 104178569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 104278569bb4SMatthew G. Knepley 104320f4b53cSBarry Smith Collective 104478569bb4SMatthew G. Knepley 104578569bb4SMatthew G. Knepley Input Parameters: 104678569bb4SMatthew G. Knepley + v - The vector to be written 104778569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 10486823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 10496823f3c5SBlaise Bourdin - offset - the location of the variable in the file 105078569bb4SMatthew G. Knepley 105120f4b53cSBarry Smith Level: beginner 105220f4b53cSBarry Smith 105378569bb4SMatthew G. Knepley Notes: 105478569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 105578569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 105678569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 105778569bb4SMatthew G. Knepley amongst all zonal variables. 105878569bb4SMatthew G. Knepley 1059c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 106078569bb4SMatthew G. Knepley */ 1061d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1062d71ae5a4SJacob Faibussowitsch { 106378569bb4SMatthew G. Knepley MPI_Comm comm; 106478569bb4SMatthew G. Knepley PetscMPIInt size; 106578569bb4SMatthew G. Knepley DM dm; 106678569bb4SMatthew G. Knepley Vec vNatural, vComp; 106722a7544eSVaclav Hapla const PetscScalar *varray; 106878569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 106978569bb4SMatthew G. Knepley PetscBool useNatural; 107078569bb4SMatthew G. Knepley IS compIS; 107178569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 107278569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 107378569bb4SMatthew G. Knepley 107478569bb4SMatthew G. Knepley PetscFunctionBegin; 10759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 10769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10779566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 10789566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 107978569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 108078569bb4SMatthew G. Knepley if (useNatural) { 1081b7352c5cSAlexis Marboeuf PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 10829566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 10839566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 108478569bb4SMatthew G. Knepley } else { 108578569bb4SMatthew G. Knepley vNatural = v; 108678569bb4SMatthew G. Knepley } 10876823f3c5SBlaise Bourdin 108878569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 108978569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 109078569bb4SMatthew G. Knepley We assume that they are stored sequentially 109178569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1092a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 109378569bb4SMatthew G. Knepley to figure out what to save where. */ 109478569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 10959566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1096792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 109778569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 109878569bb4SMatthew G. Knepley ex_block block; 109978569bb4SMatthew G. Knepley 110078569bb4SMatthew G. Knepley block.id = csID[set]; 110178569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1102792fecdfSBarry Smith PetscCallExternal(ex_get_block_param, exoid, &block); 110378569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 110478569bb4SMatthew G. Knepley } 11059566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 11069566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 11079566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 110878569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 110978569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 111078569bb4SMatthew G. Knepley 111178569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 111278569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 111378569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 111478569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 111578569bb4SMatthew G. Knepley if (bs == 1) { 11169566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 1117792fecdfSBarry Smith PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]); 11189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 111978569bb4SMatthew G. Knepley } else { 112078569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 11219566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 11229566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 11239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 1124792fecdfSBarry Smith 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)]); 11259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 11269566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 112778569bb4SMatthew G. Knepley } 112878569bb4SMatthew G. Knepley } 112978569bb4SMatthew G. Knepley csxs += csSize[set]; 113078569bb4SMatthew G. Knepley } 11319566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 11329566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 1133b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(VecDestroy(&vNatural)); 11343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113578569bb4SMatthew G. Knepley } 113678569bb4SMatthew G. Knepley 113778569bb4SMatthew G. Knepley /* 113878569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 113978569bb4SMatthew G. Knepley 114020f4b53cSBarry Smith Collective 114178569bb4SMatthew G. Knepley 114278569bb4SMatthew G. Knepley Input Parameters: 114378569bb4SMatthew G. Knepley + v - The vector to be written 114478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 11456823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 11466823f3c5SBlaise Bourdin - offset - the location of the variable in the file 114778569bb4SMatthew G. Knepley 114820f4b53cSBarry Smith Level: beginner 114920f4b53cSBarry Smith 115078569bb4SMatthew G. Knepley Notes: 115178569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 115278569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 115378569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 115478569bb4SMatthew G. Knepley amongst all zonal variables. 115578569bb4SMatthew G. Knepley 1156db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 115778569bb4SMatthew G. Knepley */ 1158d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1159d71ae5a4SJacob Faibussowitsch { 116078569bb4SMatthew G. Knepley MPI_Comm comm; 116178569bb4SMatthew G. Knepley PetscMPIInt size; 116278569bb4SMatthew G. Knepley DM dm; 116378569bb4SMatthew G. Knepley Vec vNatural, vComp; 116422a7544eSVaclav Hapla PetscScalar *varray; 116578569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 116678569bb4SMatthew G. Knepley PetscBool useNatural; 116778569bb4SMatthew G. Knepley IS compIS; 116878569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 116978569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 117078569bb4SMatthew G. Knepley 117178569bb4SMatthew G. Knepley PetscFunctionBegin; 11729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 11739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11749566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 11759566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 117678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1177b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1178ad540459SPierre Jolivet else vNatural = v; 11796823f3c5SBlaise Bourdin 118078569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 118178569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 118278569bb4SMatthew G. Knepley We assume that they are stored sequentially 118378569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1184a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 118578569bb4SMatthew G. Knepley to figure out what to save where. */ 118678569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 11879566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1188792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 118978569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 119078569bb4SMatthew G. Knepley ex_block block; 119178569bb4SMatthew G. Knepley 119278569bb4SMatthew G. Knepley block.id = csID[set]; 119378569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1194792fecdfSBarry Smith PetscCallExternal(ex_get_block_param, exoid, &block); 119578569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 119678569bb4SMatthew G. Knepley } 11979566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 11989566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 11999566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 120078569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 120178569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 120278569bb4SMatthew G. Knepley 120378569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 120478569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 120578569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 120678569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 120778569bb4SMatthew G. Knepley if (bs == 1) { 12089566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 1209792fecdfSBarry Smith PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]); 12109566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 121178569bb4SMatthew G. Knepley } else { 121278569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 12139566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 12149566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 12159566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 1216792fecdfSBarry Smith 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)]); 12179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 12189566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 121978569bb4SMatthew G. Knepley } 122078569bb4SMatthew G. Knepley } 122178569bb4SMatthew G. Knepley csxs += csSize[set]; 122278569bb4SMatthew G. Knepley } 12239566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 12249566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 122578569bb4SMatthew G. Knepley if (useNatural) { 12269566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 12279566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1228b7352c5cSAlexis Marboeuf PetscCall(VecDestroy(&vNatural)); 122978569bb4SMatthew G. Knepley } 12303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123178569bb4SMatthew G. Knepley } 1232b53c8456SSatish Balay #endif 123378569bb4SMatthew G. Knepley 12341e50132fSMatthew G. Knepley /*@ 1235a1cb98faSBarry Smith PetscViewerExodusIIGetId - Get the file id of the `PETSCVIEWEREXODUSII` file 12361e50132fSMatthew G. Knepley 123720f4b53cSBarry Smith Logically Collective 12381e50132fSMatthew G. Knepley 12391e50132fSMatthew G. Knepley Input Parameter: 1240a1cb98faSBarry Smith . viewer - the `PetscViewer` 12411e50132fSMatthew G. Knepley 12421e50132fSMatthew G. Knepley Output Parameter: 1243d8d19677SJose E. Roman . exoid - The ExodusII file id 12441e50132fSMatthew G. Knepley 12451e50132fSMatthew G. Knepley Level: intermediate 12461e50132fSMatthew G. Knepley 1247a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 12481e50132fSMatthew G. Knepley @*/ 1249d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) 1250d71ae5a4SJacob Faibussowitsch { 12511e50132fSMatthew G. Knepley PetscFunctionBegin; 12521e50132fSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1253cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetId_C", (PetscViewer, int *), (viewer, exoid)); 12543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12556823f3c5SBlaise Bourdin } 12566823f3c5SBlaise Bourdin 12576823f3c5SBlaise Bourdin /*@ 12586823f3c5SBlaise Bourdin PetscViewerExodusIISetOrder - Set the elements order in the exodusII file. 12596823f3c5SBlaise Bourdin 12606823f3c5SBlaise Bourdin Collective 12616823f3c5SBlaise Bourdin 12626823f3c5SBlaise Bourdin Input Parameters: 1263a1cb98faSBarry Smith + viewer - the `PETSCVIEWEREXODUSII` viewer 126498a6a78aSPierre Jolivet - order - elements order 12656823f3c5SBlaise Bourdin 12666823f3c5SBlaise Bourdin Output Parameter: 12676823f3c5SBlaise Bourdin 12686823f3c5SBlaise Bourdin Level: beginner 12696823f3c5SBlaise Bourdin 1270a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 12716823f3c5SBlaise Bourdin @*/ 1272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order) 1273d71ae5a4SJacob Faibussowitsch { 12746823f3c5SBlaise Bourdin PetscFunctionBegin; 12756823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1276cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerSetOrder_C", (PetscViewer, PetscInt), (viewer, order)); 12773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12786823f3c5SBlaise Bourdin } 12796823f3c5SBlaise Bourdin 12806823f3c5SBlaise Bourdin /*@ 12816823f3c5SBlaise Bourdin PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file. 12826823f3c5SBlaise Bourdin 12836823f3c5SBlaise Bourdin Collective 12846823f3c5SBlaise Bourdin 12856823f3c5SBlaise Bourdin Input Parameters: 1286a1cb98faSBarry Smith + viewer - the `PETSCVIEWEREXODUSII` viewer 128798a6a78aSPierre Jolivet - order - elements order 12886823f3c5SBlaise Bourdin 12896823f3c5SBlaise Bourdin Output Parameter: 12906823f3c5SBlaise Bourdin 12916823f3c5SBlaise Bourdin Level: beginner 12926823f3c5SBlaise Bourdin 1293a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 12946823f3c5SBlaise Bourdin @*/ 1295d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order) 1296d71ae5a4SJacob Faibussowitsch { 12976823f3c5SBlaise Bourdin PetscFunctionBegin; 12986823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1299cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetOrder_C", (PetscViewer, PetscInt *), (viewer, order)); 13003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13016823f3c5SBlaise Bourdin } 13026823f3c5SBlaise Bourdin 13036823f3c5SBlaise Bourdin /*@C 13046823f3c5SBlaise Bourdin PetscViewerExodusIIOpen - Opens a file for ExodusII input/output. 13056823f3c5SBlaise Bourdin 13066823f3c5SBlaise Bourdin Collective 13076823f3c5SBlaise Bourdin 13086823f3c5SBlaise Bourdin Input Parameters: 13096823f3c5SBlaise Bourdin + comm - MPI communicator 13106823f3c5SBlaise Bourdin . name - name of file 13116823f3c5SBlaise Bourdin - type - type of file 1312a1cb98faSBarry Smith .vb 1313a1cb98faSBarry Smith FILE_MODE_WRITE - create new file for binary output 1314a1cb98faSBarry Smith FILE_MODE_READ - open existing file for binary input 1315a1cb98faSBarry Smith FILE_MODE_APPEND - open existing file for binary output 1316a1cb98faSBarry Smith .ve 13176823f3c5SBlaise Bourdin 13186823f3c5SBlaise Bourdin Output Parameter: 1319a1cb98faSBarry Smith . exo - `PETSCVIEWEREXODUSII` `PetscViewer` for Exodus II input/output to use with the specified file 13206823f3c5SBlaise Bourdin 13216823f3c5SBlaise Bourdin Level: beginner 13226823f3c5SBlaise Bourdin 1323a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 1324db781477SPatrick Sanan `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 13256823f3c5SBlaise Bourdin @*/ 1326d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) 1327d71ae5a4SJacob Faibussowitsch { 13286823f3c5SBlaise Bourdin PetscFunctionBegin; 13299566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, exo)); 13309566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII)); 13319566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*exo, type)); 13329566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*exo, name)); 13339566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*exo)); 13343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13351e50132fSMatthew G. Knepley } 13361e50132fSMatthew G. Knepley 133733751fbdSMichael Lange /*@C 1338a1cb98faSBarry Smith DMPlexCreateExodusFromFile - Create a `DMPLEX` mesh from an ExodusII file. 133933751fbdSMichael Lange 1340d083f849SBarry Smith Collective 134133751fbdSMichael Lange 134233751fbdSMichael Lange Input Parameters: 134333751fbdSMichael Lange + comm - The MPI communicator 134433751fbdSMichael Lange . filename - The name of the ExodusII file 134533751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 134633751fbdSMichael Lange 134733751fbdSMichael Lange Output Parameter: 1348a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 134933751fbdSMichael Lange 135033751fbdSMichael Lange Level: beginner 135133751fbdSMichael Lange 1352*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()` 135333751fbdSMichael Lange @*/ 1354d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1355d71ae5a4SJacob Faibussowitsch { 135633751fbdSMichael Lange PetscMPIInt rank; 135733751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1358e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 135933751fbdSMichael Lange float version; 136033751fbdSMichael Lange #endif 136133751fbdSMichael Lange 136233751fbdSMichael Lange PetscFunctionBegin; 136333751fbdSMichael Lange PetscValidCharPointer(filename, 2); 13649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 136533751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1366dd400576SPatrick Sanan if (rank == 0) { 136733751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 136808401ef6SPierre Jolivet PetscCheck(exoid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 136933751fbdSMichael Lange } 13709566063dSJacob Faibussowitsch PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm)); 137148a46eb9SPierre Jolivet if (rank == 0) PetscCallExternal(ex_close, exoid); 13723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 137333751fbdSMichael Lange #else 137433751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 137533751fbdSMichael Lange #endif 137633751fbdSMichael Lange } 137733751fbdSMichael Lange 13788f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 1379d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct) 1380d71ae5a4SJacob Faibussowitsch { 13818f861fbcSMatthew G. Knepley PetscBool flg; 13828f861fbcSMatthew G. Knepley 13838f861fbcSMatthew G. Knepley PetscFunctionBegin; 13848f861fbcSMatthew G. Knepley *ct = DM_POLYTOPE_UNKNOWN; 13859566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI", &flg)); 13869371c9d4SSatish Balay if (flg) { 13879371c9d4SSatish Balay *ct = DM_POLYTOPE_TRIANGLE; 13889371c9d4SSatish Balay goto done; 13899371c9d4SSatish Balay } 13909566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI3", &flg)); 13919371c9d4SSatish Balay if (flg) { 13929371c9d4SSatish Balay *ct = DM_POLYTOPE_TRIANGLE; 13939371c9d4SSatish Balay goto done; 13949371c9d4SSatish Balay } 13959566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD", &flg)); 13969371c9d4SSatish Balay if (flg) { 13979371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 13989371c9d4SSatish Balay goto done; 13999371c9d4SSatish Balay } 14009566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg)); 14019371c9d4SSatish Balay if (flg) { 14029371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 14039371c9d4SSatish Balay goto done; 14049371c9d4SSatish Balay } 14059566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg)); 14069371c9d4SSatish Balay if (flg) { 14079371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 14089371c9d4SSatish Balay goto done; 14099371c9d4SSatish Balay } 14109566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TETRA", &flg)); 14119371c9d4SSatish Balay if (flg) { 14129371c9d4SSatish Balay *ct = DM_POLYTOPE_TETRAHEDRON; 14139371c9d4SSatish Balay goto done; 14149371c9d4SSatish Balay } 14159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TET4", &flg)); 14169371c9d4SSatish Balay if (flg) { 14179371c9d4SSatish Balay *ct = DM_POLYTOPE_TETRAHEDRON; 14189371c9d4SSatish Balay goto done; 14199371c9d4SSatish Balay } 14209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg)); 14219371c9d4SSatish Balay if (flg) { 14229371c9d4SSatish Balay *ct = DM_POLYTOPE_TRI_PRISM; 14239371c9d4SSatish Balay goto done; 14249371c9d4SSatish Balay } 14259566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX", &flg)); 14269371c9d4SSatish Balay if (flg) { 14279371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14289371c9d4SSatish Balay goto done; 14299371c9d4SSatish Balay } 14309566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX8", &flg)); 14319371c9d4SSatish Balay if (flg) { 14329371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14339371c9d4SSatish Balay goto done; 14349371c9d4SSatish Balay } 14359566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg)); 14369371c9d4SSatish Balay if (flg) { 14379371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14389371c9d4SSatish Balay goto done; 14399371c9d4SSatish Balay } 144098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type); 14418f861fbcSMatthew G. Knepley done: 14423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14438f861fbcSMatthew G. Knepley } 14448f861fbcSMatthew G. Knepley #endif 14458f861fbcSMatthew G. Knepley 1446552f7358SJed Brown /*@ 1447a1cb98faSBarry Smith DMPlexCreateExodus - Create a `DMPLEX` mesh from an ExodusII file ID. 1448552f7358SJed Brown 1449d083f849SBarry Smith Collective 1450552f7358SJed Brown 1451552f7358SJed Brown Input Parameters: 1452552f7358SJed Brown + comm - The MPI communicator 1453552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 1454552f7358SJed Brown - interpolate - Create faces and edges in the mesh 1455552f7358SJed Brown 1456552f7358SJed Brown Output Parameter: 1457a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 1458552f7358SJed Brown 1459552f7358SJed Brown Level: beginner 1460552f7358SJed Brown 1461*1cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMPLEX`, `DMCreate()` 1462552f7358SJed Brown @*/ 1463d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 1464d71ae5a4SJacob Faibussowitsch { 1465552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 1466552f7358SJed Brown PetscMPIInt num_proc, rank; 1467ae9eebabSLisandro Dalcin DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL; 1468552f7358SJed Brown PetscSection coordSection; 1469552f7358SJed Brown Vec coordinates; 1470552f7358SJed Brown PetscScalar *coords; 1471552f7358SJed Brown PetscInt coordSize, v; 1472552f7358SJed Brown /* Read from ex_get_init() */ 1473552f7358SJed Brown char title[PETSC_MAX_PATH_LEN + 1]; 14745f80ce2aSJacob Faibussowitsch int dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0; 1475552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 1476552f7358SJed Brown #endif 1477552f7358SJed Brown 1478552f7358SJed Brown PetscFunctionBegin; 1479552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 14809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 14819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 14829566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 14839566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1484a5b23f4aSJose E. Roman /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */ 1485dd400576SPatrick Sanan if (rank == 0) { 14869566063dSJacob Faibussowitsch PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1)); 1487792fecdfSBarry Smith PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 148828b400f6SJacob Faibussowitsch PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set"); 1489552f7358SJed Brown } 14909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm)); 14919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 14929566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, title)); 14939566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 1494412e9a14SMatthew G. Knepley /* We do not want this label automatically computed, instead we compute it here */ 14959566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "celltype")); 1496552f7358SJed Brown 1497552f7358SJed Brown /* Read cell sets information */ 1498dd400576SPatrick Sanan if (rank == 0) { 1499552f7358SJed Brown PetscInt *cone; 1500e8f6893fSMatthew G. Knepley int c, cs, ncs, c_loc, v, v_loc; 1501552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 1502e8f6893fSMatthew G. Knepley int *cs_id, *cs_order; 1503552f7358SJed Brown /* Read from ex_get_elem_block() */ 1504552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN + 1]; 1505e8f6893fSMatthew G. Knepley int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr; 1506552f7358SJed Brown /* Read from ex_get_elem_conn() */ 1507552f7358SJed Brown int *cs_connect; 1508552f7358SJed Brown 1509552f7358SJed Brown /* Get cell sets IDs */ 15109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order)); 1511792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id); 1512552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 1513552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 1514e8f6893fSMatthew G. Knepley /* Check for a hybrid mesh */ 1515e8f6893fSMatthew G. Knepley for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) { 15168f861fbcSMatthew G. Knepley DMPolytopeType ct; 15178f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 15188f861fbcSMatthew G. Knepley 15199566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1520792fecdfSBarry Smith PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 15219566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 15228f861fbcSMatthew G. Knepley dim = PetscMax(dim, DMPolytopeTypeGetDim(ct)); 1523792fecdfSBarry Smith PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 15248f861fbcSMatthew G. Knepley switch (ct) { 15258f861fbcSMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 1526e8f6893fSMatthew G. Knepley cs_order[cs] = cs; 1527e8f6893fSMatthew G. Knepley ++num_hybrid; 1528e8f6893fSMatthew G. Knepley break; 1529e8f6893fSMatthew G. Knepley default: 1530e8f6893fSMatthew G. Knepley for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1]; 1531e8f6893fSMatthew G. Knepley cs_order[cs - num_hybrid] = cs; 1532e8f6893fSMatthew G. Knepley } 1533e8f6893fSMatthew G. Knepley } 1534552f7358SJed Brown /* First set sizes */ 1535e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 15368f861fbcSMatthew G. Knepley DMPolytopeType ct; 15378f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 1538e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 15398f861fbcSMatthew G. Knepley 15409566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1541792fecdfSBarry Smith PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 15429566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 1543792fecdfSBarry Smith PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 1544552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 15459566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell)); 15469566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, c, ct)); 1547552f7358SJed Brown } 1548552f7358SJed Brown } 15499566063dSJacob Faibussowitsch for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 15509566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 1551e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1552e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 1553792fecdfSBarry Smith PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 15549566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone)); 1555792fecdfSBarry Smith PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL); 1556eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 1557552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1558636e64ffSStefano Zampini DMPolytopeType ct; 1559636e64ffSStefano Zampini 1560ad540459SPierre Jolivet for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1; 15619566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(*dm, c, &ct)); 15629566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCell(ct, cone)); 15639566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, cone)); 15649566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs])); 1565552f7358SJed Brown } 15669566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_connect, cone)); 1567552f7358SJed Brown } 15689566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_id, cs_order)); 1569552f7358SJed Brown } 15708f861fbcSMatthew G. Knepley { 15718f861fbcSMatthew G. Knepley PetscInt ints[] = {dim, dimEmbed}; 15728f861fbcSMatthew G. Knepley 15739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm)); 15749566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, ints[0])); 15759566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, ints[1])); 15768f861fbcSMatthew G. Knepley dim = ints[0]; 15778f861fbcSMatthew G. Knepley dimEmbed = ints[1]; 15788f861fbcSMatthew G. Knepley } 15799566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 15809566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1581552f7358SJed Brown if (interpolate) { 15825fd9971aSMatthew G. Knepley DM idm; 1583552f7358SJed Brown 15849566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 15859566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1586552f7358SJed Brown *dm = idm; 1587552f7358SJed Brown } 1588552f7358SJed Brown 1589552f7358SJed Brown /* Create vertex set label */ 1590dd400576SPatrick Sanan if (rank == 0 && (num_vs > 0)) { 1591552f7358SJed Brown int vs, v; 1592552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 1593552f7358SJed Brown int *vs_id; 1594552f7358SJed Brown /* Read from ex_get_node_set_param() */ 1595f958083aSBlaise Bourdin int num_vertex_in_set; 1596552f7358SJed Brown /* Read from ex_get_node_set() */ 1597552f7358SJed Brown int *vs_vertex_list; 1598552f7358SJed Brown 1599552f7358SJed Brown /* Get vertex set ids */ 16009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vs, &vs_id)); 1601792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id); 1602552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 1603792fecdfSBarry Smith PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 16049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list)); 1605792fecdfSBarry Smith PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL); 160648a46eb9SPierre Jolivet 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])); 16079566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_vertex_list)); 1608552f7358SJed Brown } 16099566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_id)); 1610552f7358SJed Brown } 1611552f7358SJed Brown /* Read coordinates */ 16129566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 16139566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 16149566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed)); 16159566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 1616552f7358SJed Brown for (v = numCells; v < numCells + numVertices; ++v) { 16179566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed)); 16189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed)); 1619552f7358SJed Brown } 16209566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 16219566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 16229566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 16239566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 16249566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 16259566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, dimEmbed)); 16269566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 16279566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 1628dd400576SPatrick Sanan if (rank == 0) { 1629e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 1630552f7358SJed Brown 16319566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z)); 1632792fecdfSBarry Smith PetscCallExternal(ex_get_coord, exoid, x, y, z); 16338f861fbcSMatthew G. Knepley if (dimEmbed > 0) { 16348f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v]; 16350d644c17SKarl Rupp } 16368f861fbcSMatthew G. Knepley if (dimEmbed > 1) { 16378f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v]; 16380d644c17SKarl Rupp } 16398f861fbcSMatthew G. Knepley if (dimEmbed > 2) { 16408f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v]; 16410d644c17SKarl Rupp } 16429566063dSJacob Faibussowitsch PetscCall(PetscFree3(x, y, z)); 1643552f7358SJed Brown } 16449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 16459566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 16469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1647552f7358SJed Brown 1648552f7358SJed Brown /* Create side set label */ 1649dd400576SPatrick Sanan if (rank == 0 && interpolate && (num_fs > 0)) { 1650552f7358SJed Brown int fs, f, voff; 1651552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 1652552f7358SJed Brown int *fs_id; 1653552f7358SJed Brown /* Read from ex_get_side_set_param() */ 1654f958083aSBlaise Bourdin int num_side_in_set; 1655552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 1656552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 1657ef073283Smichael_afanasiev /* Read side set labels */ 1658ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH + 1]; 1659034d25ccSMatthew G. Knepley size_t fs_name_len; 1660552f7358SJed Brown 1661552f7358SJed Brown /* Get side set ids */ 16629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_fs, &fs_id)); 1663792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id); 1664552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 1665792fecdfSBarry Smith PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL); 16669566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list)); 1667792fecdfSBarry Smith PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list); 1668ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 1669ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1670034d25ccSMatthew G. Knepley if (!fs_name_err) { 1671034d25ccSMatthew G. Knepley PetscCall(PetscStrlen(fs_name, &fs_name_len)); 1672034d25ccSMatthew G. Knepley if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH)); 1673034d25ccSMatthew G. Knepley } 1674552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 16750298fd71SBarry Smith const PetscInt *faces = NULL; 1676552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1677552f7358SJed Brown PetscInt faceVertices[4], v; 1678552f7358SJed Brown 1679197f6770SJed Brown PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize); 1680ad540459SPierre Jolivet for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1; 16819566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1682197f6770SJed Brown PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces); 16839566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs])); 1684ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 1685034d25ccSMatthew G. Knepley if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs])); 16869566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1687552f7358SJed Brown } 16889566063dSJacob Faibussowitsch PetscCall(PetscFree2(fs_vertex_count_list, fs_vertex_list)); 1689552f7358SJed Brown } 16909566063dSJacob Faibussowitsch PetscCall(PetscFree(fs_id)); 1691552f7358SJed Brown } 1692ae9eebabSLisandro Dalcin 1693ae9eebabSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 16949371c9d4SSatish Balay enum { 16959371c9d4SSatish Balay n = 3 16969371c9d4SSatish Balay }; 1697ae9eebabSLisandro Dalcin PetscBool flag[n]; 1698ae9eebabSLisandro Dalcin 1699ae9eebabSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 1700ae9eebabSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 1701ae9eebabSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 17029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 17039566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 17049566063dSJacob Faibussowitsch if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 17059566063dSJacob Faibussowitsch if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 1706ae9eebabSLisandro Dalcin } 17073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1708552f7358SJed Brown #else 1709552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1710552f7358SJed Brown #endif 1711552f7358SJed Brown } 1712