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 28760225df5SJacob Faibussowitsch - 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 { 310207ab81fSDavid Andrs SEGMENT, 3119371c9d4SSatish Balay TRI, 3129371c9d4SSatish Balay QUAD, 3139371c9d4SSatish Balay TET, 3149371c9d4SSatish Balay HEX 3159371c9d4SSatish Balay }; 31678569bb4SMatthew G. Knepley MPI_Comm comm; 3176823f3c5SBlaise Bourdin PetscInt degree; /* the order of the mesh */ 31878569bb4SMatthew G. Knepley /* Connectivity Variables */ 31978569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 32078569bb4SMatthew G. Knepley /* Cell Sets */ 32178569bb4SMatthew G. Knepley DMLabel csLabel; 32278569bb4SMatthew G. Knepley IS csIS; 32378569bb4SMatthew G. Knepley const PetscInt *csIdx; 32478569bb4SMatthew G. Knepley PetscInt num_cs, cs; 32578569bb4SMatthew G. Knepley enum ElemType *type; 326fe209ef9SBlaise Bourdin PetscBool hasLabel; 32778569bb4SMatthew G. Knepley /* Coordinate Variables */ 32878569bb4SMatthew G. Knepley DM cdm; 3296823f3c5SBlaise Bourdin PetscSection coordSection; 33078569bb4SMatthew G. Knepley Vec coord; 33178569bb4SMatthew G. Knepley PetscInt **nodes; 332eae8a387SMatthew G. Knepley PetscInt depth, d, dim, skipCells = 0; 33378569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 33478569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 33578569bb4SMatthew G. Knepley PetscMPIInt rank, size; 33678569bb4SMatthew G. Knepley const char *dmName; 337207ab81fSDavid Andrs PetscInt nodesLineP1[4] = {2, 0, 0, 0}; 338207ab81fSDavid Andrs PetscInt nodesLineP2[4] = {2, 0, 0, 1}; 339fe209ef9SBlaise Bourdin PetscInt nodesTriP1[4] = {3, 0, 0, 0}; 340fe209ef9SBlaise Bourdin PetscInt nodesTriP2[4] = {3, 3, 0, 0}; 341fe209ef9SBlaise Bourdin PetscInt nodesQuadP1[4] = {4, 0, 0, 0}; 342fe209ef9SBlaise Bourdin PetscInt nodesQuadP2[4] = {4, 4, 0, 1}; 343fe209ef9SBlaise Bourdin PetscInt nodesTetP1[4] = {4, 0, 0, 0}; 344fe209ef9SBlaise Bourdin PetscInt nodesTetP2[4] = {4, 6, 0, 0}; 345fe209ef9SBlaise Bourdin PetscInt nodesHexP1[4] = {8, 0, 0, 0}; 346fe209ef9SBlaise Bourdin PetscInt nodesHexP2[4] = {8, 12, 6, 1}; 3476823f3c5SBlaise Bourdin int CPU_word_size, IO_word_size, EXO_mode; 3486823f3c5SBlaise Bourdin float EXO_version; 3496823f3c5SBlaise Bourdin 3506823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data; 3516823f3c5SBlaise Bourdin 35278569bb4SMatthew G. Knepley PetscFunctionBegin; 3539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3566823f3c5SBlaise Bourdin 3576823f3c5SBlaise Bourdin /* 3586823f3c5SBlaise Bourdin Creating coordSection is a collective operation so we do it somewhat out of sequence 3596823f3c5SBlaise Bourdin */ 3609566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &coordSection)); 3619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(dm)); 362bdbfaa5dSBlaise Bourdin /* 363bdbfaa5dSBlaise Bourdin Check that all points are on rank 0 since we don't know how to save distributed DM in exodus format 364bdbfaa5dSBlaise Bourdin */ 365bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 366bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 367bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 368bdbfaa5dSBlaise Bourdin PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 369bdbfaa5dSBlaise Bourdin numCells = cEnd - cStart; 370bdbfaa5dSBlaise Bourdin numEdges = eEnd - eStart; 371bdbfaa5dSBlaise Bourdin numVertices = vEnd - vStart; 372bdbfaa5dSBlaise Bourdin PetscCheck(!(rank && (numCells || numEdges || numVertices)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Writing distributed DM in exodusII format not supported"); 373c5853193SPierre Jolivet if (rank == 0) { 3746823f3c5SBlaise Bourdin switch (exo->btype) { 3756823f3c5SBlaise Bourdin case FILE_MODE_READ: 3766823f3c5SBlaise Bourdin case FILE_MODE_APPEND: 3776823f3c5SBlaise Bourdin case FILE_MODE_UPDATE: 3786823f3c5SBlaise Bourdin case FILE_MODE_APPEND_UPDATE: 3796823f3c5SBlaise Bourdin /* exodusII does not allow writing geometry to an existing file */ 38098921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename); 3816823f3c5SBlaise Bourdin case FILE_MODE_WRITE: 3826823f3c5SBlaise Bourdin /* Create an empty file if one already exists*/ 3836823f3c5SBlaise Bourdin EXO_mode = EX_CLOBBER; 3846823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 3856823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 3866823f3c5SBlaise Bourdin #endif 3876823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 3886823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 3896823f3c5SBlaise Bourdin exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size); 39008401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename); 3916823f3c5SBlaise Bourdin 3926823f3c5SBlaise Bourdin break; 393d71ae5a4SJacob Faibussowitsch default: 394d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3956823f3c5SBlaise Bourdin } 3966823f3c5SBlaise Bourdin 39778569bb4SMatthew G. Knepley /* --- Get DM info --- */ 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm, &dmName)); 3999566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 4009566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 4019566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 4029371c9d4SSatish Balay if (depth == 3) { 4039371c9d4SSatish Balay numFaces = fEnd - fStart; 4049371c9d4SSatish Balay } else { 4059371c9d4SSatish Balay numFaces = 0; 4069371c9d4SSatish Balay } 4079566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs)); 4089566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs)); 4099566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs)); 4109566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coord)); 4119566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 41278569bb4SMatthew G. Knepley if (num_cs > 0) { 4139566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel)); 4149566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(csLabel, &csIS)); 4159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(csIS, &csIdx)); 41678569bb4SMatthew G. Knepley } 4179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &nodes)); 41878569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 4199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &type)); 42078569bb4SMatthew G. Knepley numNodes = numVertices; 4216823f3c5SBlaise Bourdin 4229566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetOrder(viewer, °ree)); 423ad540459SPierre Jolivet if (degree == 2) numNodes += numEdges; 42478569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 42578569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 42678569bb4SMatthew G. Knepley IS stratumIS; 42778569bb4SMatthew G. Knepley const PetscInt *cells; 42878569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 42978569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 43078569bb4SMatthew G. Knepley 4319566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 4329566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 4339566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 4349566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 43578569bb4SMatthew G. Knepley switch (dim) { 436207ab81fSDavid Andrs case 1: 437207ab81fSDavid Andrs if (closureSize == 2 * dim) { 438207ab81fSDavid Andrs type[cs] = SEGMENT; 439207ab81fSDavid Andrs } 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); 44078569bb4SMatthew G. Knepley case 2: 4419371c9d4SSatish Balay if (closureSize == 3 * dim) { 4429371c9d4SSatish Balay type[cs] = TRI; 4439371c9d4SSatish Balay } else if (closureSize == 4 * dim) { 4449371c9d4SSatish Balay type[cs] = QUAD; 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; 44778569bb4SMatthew G. Knepley case 3: 4489371c9d4SSatish Balay if (closureSize == 4 * dim) { 4499371c9d4SSatish Balay type[cs] = TET; 4509371c9d4SSatish Balay } else if (closureSize == 8 * dim) { 4519371c9d4SSatish Balay type[cs] = HEX; 4529371c9d4SSatish 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); 45378569bb4SMatthew G. Knepley break; 454d71ae5a4SJacob Faibussowitsch default: 455d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim); 45678569bb4SMatthew G. Knepley } 457207ab81fSDavid Andrs if ((degree == 2) && (type[cs] == SEGMENT)) numNodes += csSize; 458ad540459SPierre Jolivet if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize; 4599371c9d4SSatish Balay if ((degree == 2) && (type[cs] == HEX)) { 4609371c9d4SSatish Balay numNodes += csSize; 4619371c9d4SSatish Balay numNodes += numFaces; 4629371c9d4SSatish Balay } 4639566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 46478569bb4SMatthew G. Knepley /* Set nodes and Element type */ 465207ab81fSDavid Andrs if (type[cs] == SEGMENT) { 466207ab81fSDavid Andrs if (degree == 1) nodes[cs] = nodesLineP1; 467207ab81fSDavid Andrs else if (degree == 2) nodes[cs] = nodesLineP2; 468207ab81fSDavid Andrs } else if (type[cs] == TRI) { 46978569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 47078569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 47178569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 47278569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 47378569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 47478569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 47578569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 47678569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 47778569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 47878569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 47978569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 48078569bb4SMatthew G. Knepley } 48178569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 48278569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3] * csSize; 48378569bb4SMatthew G. Knepley 4849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 4859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 48678569bb4SMatthew G. Knepley } 487bdbfaa5dSBlaise Bourdin if (num_cs) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs); 48878569bb4SMatthew G. Knepley /* --- Connectivity --- */ 48978569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 49078569bb4SMatthew G. Knepley IS stratumIS; 49178569bb4SMatthew G. Knepley const PetscInt *cells; 492fe209ef9SBlaise Bourdin PetscInt *connect, off = 0; 493eae8a387SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 49478569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 49578569bb4SMatthew G. Knepley char *elem_type = NULL; 496207ab81fSDavid Andrs char elem_type_bar2[] = "BAR2", elem_type_bar3[] = "BAR3"; 49778569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 49878569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 49978569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 50078569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 50178569bb4SMatthew G. Knepley 5029566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 5039566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 5049566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 50578569bb4SMatthew G. Knepley /* Set Element type */ 506207ab81fSDavid Andrs if (type[cs] == SEGMENT) { 507207ab81fSDavid Andrs if (degree == 1) elem_type = elem_type_bar2; 508207ab81fSDavid Andrs else if (degree == 2) elem_type = elem_type_bar3; 509207ab81fSDavid Andrs } else if (type[cs] == TRI) { 51078569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 51178569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 51278569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 51378569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 51478569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 51578569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 51678569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 51778569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 51878569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 51978569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 52078569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 52178569bb4SMatthew G. Knepley } 52278569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect)); 524792fecdfSBarry Smith PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1); 52578569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 52678569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 52778569bb4SMatthew G. Knepley if (depth > 1) { 52878569bb4SMatthew G. Knepley if (dim == 2) { 5299566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure)); 53078569bb4SMatthew G. Knepley } else if (dim == 3) { 53178569bb4SMatthew G. Knepley PetscInt *closure = NULL; 53278569bb4SMatthew G. Knepley 5339566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure)); 5349566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 53578569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 5369566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 53778569bb4SMatthew G. Knepley } 53878569bb4SMatthew G. Knepley } 53978569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 54078569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 54178569bb4SMatthew G. Knepley PetscInt *closure = NULL; 54278569bb4SMatthew G. Knepley PetscInt temp, i; 54378569bb4SMatthew G. Knepley 5449566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 54578569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 54678569bb4SMatthew G. Knepley if (i < nodes[cs][0]) { /* Vertices */ 547fe209ef9SBlaise Bourdin connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1; 548fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 54978569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */ 550fe209ef9SBlaise Bourdin connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1; 551fe209ef9SBlaise Bourdin if (nodes[cs][2] == 0) connect[i + off] -= numFaces; 552fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 55378569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */ 554fe209ef9SBlaise Bourdin connect[i + off] = closure[0] + 1; 555fe209ef9SBlaise Bourdin connect[i + off] -= skipCells; 55678569bb4SMatthew G. Knepley } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */ 557fe209ef9SBlaise Bourdin connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1; 558fe209ef9SBlaise Bourdin connect[i + off] -= cellsNotInConnectivity; 55978569bb4SMatthew G. Knepley } else { 560fe209ef9SBlaise Bourdin connect[i + off] = -1; 56178569bb4SMatthew G. Knepley } 56278569bb4SMatthew G. Knepley } 56378569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 56478569bb4SMatthew G. Knepley if (type[cs] == TET) { 5659371c9d4SSatish Balay temp = connect[0 + off]; 5669371c9d4SSatish Balay connect[0 + off] = connect[1 + off]; 5679371c9d4SSatish Balay connect[1 + off] = temp; 56878569bb4SMatthew G. Knepley if (degree == 2) { 5699371c9d4SSatish Balay temp = connect[5 + off]; 5709371c9d4SSatish Balay connect[5 + off] = connect[6 + off]; 5719371c9d4SSatish Balay connect[6 + off] = temp; 5729371c9d4SSatish Balay temp = connect[7 + off]; 5739371c9d4SSatish Balay connect[7 + off] = connect[8 + off]; 5749371c9d4SSatish Balay connect[8 + off] = temp; 57578569bb4SMatthew G. Knepley } 57678569bb4SMatthew G. Knepley } 57778569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 57878569bb4SMatthew G. Knepley if (type[cs] == HEX) { 5799371c9d4SSatish Balay temp = connect[1 + off]; 5809371c9d4SSatish Balay connect[1 + off] = connect[3 + off]; 5819371c9d4SSatish Balay connect[3 + off] = temp; 58278569bb4SMatthew G. Knepley if (degree == 2) { 5839371c9d4SSatish Balay temp = connect[8 + off]; 5849371c9d4SSatish Balay connect[8 + off] = connect[11 + off]; 5859371c9d4SSatish Balay connect[11 + off] = temp; 5869371c9d4SSatish Balay temp = connect[9 + off]; 5879371c9d4SSatish Balay connect[9 + off] = connect[10 + off]; 5889371c9d4SSatish Balay connect[10 + off] = temp; 5899371c9d4SSatish Balay temp = connect[16 + off]; 5909371c9d4SSatish Balay connect[16 + off] = connect[17 + off]; 5919371c9d4SSatish Balay connect[17 + off] = temp; 5929371c9d4SSatish Balay temp = connect[18 + off]; 5939371c9d4SSatish Balay connect[18 + off] = connect[19 + off]; 5949371c9d4SSatish Balay connect[19 + off] = temp; 59578569bb4SMatthew G. Knepley 5969371c9d4SSatish Balay temp = connect[12 + off]; 5979371c9d4SSatish Balay connect[12 + off] = connect[16 + off]; 5989371c9d4SSatish Balay connect[16 + off] = temp; 5999371c9d4SSatish Balay temp = connect[13 + off]; 6009371c9d4SSatish Balay connect[13 + off] = connect[17 + off]; 6019371c9d4SSatish Balay connect[17 + off] = temp; 6029371c9d4SSatish Balay temp = connect[14 + off]; 6039371c9d4SSatish Balay connect[14 + off] = connect[18 + off]; 6049371c9d4SSatish Balay connect[18 + off] = temp; 6059371c9d4SSatish Balay temp = connect[15 + off]; 6069371c9d4SSatish Balay connect[15 + off] = connect[19 + off]; 6079371c9d4SSatish Balay connect[19 + off] = temp; 60878569bb4SMatthew G. Knepley 6099371c9d4SSatish Balay temp = connect[23 + off]; 6109371c9d4SSatish Balay connect[23 + off] = connect[26 + off]; 6119371c9d4SSatish Balay connect[26 + off] = temp; 6129371c9d4SSatish Balay temp = connect[24 + off]; 6139371c9d4SSatish Balay connect[24 + off] = connect[25 + off]; 6149371c9d4SSatish Balay connect[25 + off] = temp; 6159371c9d4SSatish Balay temp = connect[25 + off]; 6169371c9d4SSatish Balay connect[25 + off] = connect[26 + off]; 6179371c9d4SSatish Balay connect[26 + off] = temp; 61878569bb4SMatthew G. Knepley } 61978569bb4SMatthew G. Knepley } 620fe209ef9SBlaise Bourdin off += connectSize; 6219566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 62278569bb4SMatthew G. Knepley } 623792fecdfSBarry Smith PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0); 62478569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0) * csSize; 6259566063dSJacob Faibussowitsch PetscCall(PetscFree(connect)); 6269566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 6279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 62878569bb4SMatthew G. Knepley } 6299566063dSJacob Faibussowitsch PetscCall(PetscFree(type)); 63078569bb4SMatthew G. Knepley /* --- Coordinates --- */ 6319566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd)); 6321e50132fSMatthew G. Knepley if (num_cs) { 63378569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 6349566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 63548a46eb9SPierre Jolivet for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0)); 63678569bb4SMatthew G. Knepley } 6371e50132fSMatthew G. Knepley } 63878569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 63978569bb4SMatthew G. Knepley IS stratumIS; 64078569bb4SMatthew G. Knepley const PetscInt *cells; 64178569bb4SMatthew G. Knepley PetscInt csSize, c; 64278569bb4SMatthew G. Knepley 6439566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 6449566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 6459566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 64648a46eb9SPierre Jolivet for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0)); 6479566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 6489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 64978569bb4SMatthew G. Knepley } 650bdbfaa5dSBlaise Bourdin if (num_cs) { 6519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(csIS, &csIdx)); 6529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS)); 65378569bb4SMatthew G. Knepley } 6549566063dSJacob Faibussowitsch PetscCall(PetscFree(nodes)); 6559566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 656bdbfaa5dSBlaise Bourdin if (numNodes) { 65778569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 658233c95e0SFrancesco Ballarin PetscScalar *closure, *cval; 659233c95e0SFrancesco Ballarin PetscReal *coords; 66078569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 66178569bb4SMatthew G. Knepley 6626823f3c5SBlaise Bourdin /* There can't be more than 24 values in the closure of a point for the coord coordSection */ 6639566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure)); 6649566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord)); 6659566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 66678569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 6679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, p, &hasDof)); 66878569bb4SMatthew G. Knepley if (hasDof) { 66978569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 67078569bb4SMatthew G. Knepley 6719566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure)); 67278569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 67378569bb4SMatthew G. Knepley cval[d] = 0.0; 67478569bb4SMatthew G. Knepley for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d]; 675233c95e0SFrancesco Ballarin coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize; 67678569bb4SMatthew G. Knepley } 67778569bb4SMatthew G. Knepley ++n; 67878569bb4SMatthew G. Knepley } 67978569bb4SMatthew G. Knepley } 680792fecdfSBarry Smith PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]); 6819566063dSJacob Faibussowitsch PetscCall(PetscFree3(coords, cval, closure)); 682792fecdfSBarry Smith PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames); 68378569bb4SMatthew G. Knepley } 6846823f3c5SBlaise Bourdin 685fe209ef9SBlaise Bourdin /* --- Node Sets/Vertex Sets --- */ 6863ba16761SJacob Faibussowitsch PetscCall(DMHasLabel(dm, "Vertex Sets", &hasLabel)); 687fe209ef9SBlaise Bourdin if (hasLabel) { 688fe209ef9SBlaise Bourdin PetscInt i, vs, vsSize; 689fe209ef9SBlaise Bourdin const PetscInt *vsIdx, *vertices; 690fe209ef9SBlaise Bourdin PetscInt *nodeList; 691fe209ef9SBlaise Bourdin IS vsIS, stratumIS; 692fe209ef9SBlaise Bourdin DMLabel vsLabel; 6939566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel)); 6949566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(vsLabel, &vsIS)); 6959566063dSJacob Faibussowitsch PetscCall(ISGetIndices(vsIS, &vsIdx)); 696fe209ef9SBlaise Bourdin for (vs = 0; vs < num_vs; ++vs) { 6979566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS)); 6989566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &vertices)); 6999566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &vsSize)); 7009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(vsSize, &nodeList)); 701ad540459SPierre Jolivet for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1; 702792fecdfSBarry Smith PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0); 703792fecdfSBarry Smith PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL); 7049566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &vertices)); 7059566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 7069566063dSJacob Faibussowitsch PetscCall(PetscFree(nodeList)); 707fe209ef9SBlaise Bourdin } 7089566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(vsIS, &vsIdx)); 7099566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vsIS)); 710fe209ef9SBlaise Bourdin } 711fe209ef9SBlaise Bourdin /* --- Side Sets/Face Sets --- */ 7129566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel)); 713fe209ef9SBlaise Bourdin if (hasLabel) { 714fe209ef9SBlaise Bourdin PetscInt i, j, fs, fsSize; 715fe209ef9SBlaise Bourdin const PetscInt *fsIdx, *faces; 716fe209ef9SBlaise Bourdin IS fsIS, stratumIS; 717fe209ef9SBlaise Bourdin DMLabel fsLabel; 718fe209ef9SBlaise Bourdin PetscInt numPoints, *points; 719fe209ef9SBlaise Bourdin PetscInt elem_list_size = 0; 720fe209ef9SBlaise Bourdin PetscInt *elem_list, *elem_ind, *side_list; 721fe209ef9SBlaise Bourdin 7229566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel)); 723fe209ef9SBlaise Bourdin /* Compute size of Node List and Element List */ 7249566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(fsLabel, &fsIS)); 7259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fsIS, &fsIdx)); 726fe209ef9SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) { 7279566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 7289566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 729fe209ef9SBlaise Bourdin elem_list_size += fsSize; 7309566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 731fe209ef9SBlaise Bourdin } 732bdbfaa5dSBlaise Bourdin if (num_fs) { 733d0609cedSBarry Smith PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list)); 734fe209ef9SBlaise Bourdin elem_ind[0] = 0; 735fe209ef9SBlaise Bourdin for (fs = 0; fs < num_fs; ++fs) { 7369566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 7379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &faces)); 7389566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 739fe209ef9SBlaise Bourdin /* Set Parameters */ 740792fecdfSBarry Smith PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0); 741fe209ef9SBlaise Bourdin /* Indices */ 742ad540459SPierre Jolivet if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize; 743fe209ef9SBlaise Bourdin 744fe209ef9SBlaise Bourdin for (i = 0; i < fsSize; ++i) { 745fe209ef9SBlaise Bourdin /* Element List */ 746fe209ef9SBlaise Bourdin points = NULL; 7479566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 748fe209ef9SBlaise Bourdin elem_list[elem_ind[fs] + i] = points[2] + 1; 7499566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 750fe209ef9SBlaise Bourdin 751fe209ef9SBlaise Bourdin /* Side List */ 752fe209ef9SBlaise Bourdin points = NULL; 7539566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 754fe209ef9SBlaise Bourdin for (j = 1; j < numPoints; ++j) { 755ad540459SPierre Jolivet if (points[j * 2] == faces[i]) break; 756fe209ef9SBlaise Bourdin } 757fe209ef9SBlaise Bourdin /* Convert HEX sides */ 758fe209ef9SBlaise Bourdin if (numPoints == 27) { 7599371c9d4SSatish Balay if (j == 1) { 7609371c9d4SSatish Balay j = 5; 7619371c9d4SSatish Balay } else if (j == 2) { 7629371c9d4SSatish Balay j = 6; 7639371c9d4SSatish Balay } else if (j == 3) { 7649371c9d4SSatish Balay j = 1; 7659371c9d4SSatish Balay } else if (j == 4) { 7669371c9d4SSatish Balay j = 3; 7679371c9d4SSatish Balay } else if (j == 5) { 7689371c9d4SSatish Balay j = 2; 7699371c9d4SSatish Balay } else if (j == 6) { 7709371c9d4SSatish Balay j = 4; 7719371c9d4SSatish Balay } 772fe209ef9SBlaise Bourdin } 773fe209ef9SBlaise Bourdin /* Convert TET sides */ 774fe209ef9SBlaise Bourdin if (numPoints == 15) { 775fe209ef9SBlaise Bourdin --j; 776ad540459SPierre Jolivet if (j == 0) j = 4; 777fe209ef9SBlaise Bourdin } 778fe209ef9SBlaise Bourdin side_list[elem_ind[fs] + i] = j; 7799566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points)); 780fe209ef9SBlaise Bourdin } 7819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &faces)); 7829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 783fe209ef9SBlaise Bourdin } 7849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fsIS, &fsIdx)); 7859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fsIS)); 786fe209ef9SBlaise Bourdin 787fe209ef9SBlaise Bourdin /* Put side sets */ 78848a46eb9SPierre 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]]); 7899566063dSJacob Faibussowitsch PetscCall(PetscFree3(elem_ind, elem_list, side_list)); 790fe209ef9SBlaise Bourdin } 791bdbfaa5dSBlaise Bourdin } 7926823f3c5SBlaise Bourdin /* 7936823f3c5SBlaise Bourdin close the exodus file 7946823f3c5SBlaise Bourdin */ 7956823f3c5SBlaise Bourdin ex_close(exo->exoid); 7966823f3c5SBlaise Bourdin exo->exoid = -1; 7976823f3c5SBlaise Bourdin } 7989566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&coordSection)); 7996823f3c5SBlaise Bourdin 8006823f3c5SBlaise Bourdin /* 8016823f3c5SBlaise Bourdin reopen the file in parallel 8026823f3c5SBlaise Bourdin */ 8036823f3c5SBlaise Bourdin EXO_mode = EX_WRITE; 8046823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 8056823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 8066823f3c5SBlaise Bourdin #endif 8076823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 8086823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 8096823f3c5SBlaise Bourdin exo->exoid = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL); 81008401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename); 8113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8126823f3c5SBlaise Bourdin } 8136823f3c5SBlaise Bourdin 8146823f3c5SBlaise Bourdin /* 8156823f3c5SBlaise Bourdin VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 8166823f3c5SBlaise Bourdin 81720f4b53cSBarry Smith Collective 8186823f3c5SBlaise Bourdin 8196823f3c5SBlaise Bourdin Input Parameters: 8206823f3c5SBlaise Bourdin + v - The vector to be written 8216823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 8226823f3c5SBlaise Bourdin 82320f4b53cSBarry Smith Level: beginner 82420f4b53cSBarry Smith 8256823f3c5SBlaise Bourdin Notes: 8266823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 8276823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 8286823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 8296823f3c5SBlaise Bourdin amongst all variables. 8306823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 8316823f3c5SBlaise Bourdin possibly corrupting the file 8326823f3c5SBlaise Bourdin 833c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 8346823f3c5SBlaise Bourdin @*/ 835d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer) 836d71ae5a4SJacob Faibussowitsch { 8376823f3c5SBlaise Bourdin DM dm; 8386823f3c5SBlaise Bourdin MPI_Comm comm; 8396823f3c5SBlaise Bourdin PetscMPIInt rank; 8406823f3c5SBlaise Bourdin 8416823f3c5SBlaise Bourdin int exoid, offsetN = 0, offsetZ = 0; 8426823f3c5SBlaise Bourdin const char *vecname; 8436823f3c5SBlaise Bourdin PetscInt step; 8446823f3c5SBlaise Bourdin 8456823f3c5SBlaise Bourdin PetscFunctionBegin; 8469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 8479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8489566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 8499566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8509566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 8516823f3c5SBlaise Bourdin 8529566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 8539566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 8549566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 8551dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 8566823f3c5SBlaise Bourdin if (offsetN > 0) { 8579566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 8586823f3c5SBlaise Bourdin } else if (offsetZ > 0) { 8599566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 86098921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 8613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8626823f3c5SBlaise Bourdin } 8636823f3c5SBlaise Bourdin 8646823f3c5SBlaise Bourdin /* 8656823f3c5SBlaise Bourdin VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 8666823f3c5SBlaise Bourdin 86720f4b53cSBarry Smith Collective 8686823f3c5SBlaise Bourdin 8696823f3c5SBlaise Bourdin Input Parameters: 8706823f3c5SBlaise Bourdin + v - The vector to be written 8716823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 8726823f3c5SBlaise Bourdin 87320f4b53cSBarry Smith Level: beginner 87420f4b53cSBarry Smith 8756823f3c5SBlaise Bourdin Notes: 8766823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 8776823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 8786823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 8796823f3c5SBlaise Bourdin amongst all variables. 8806823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 8816823f3c5SBlaise Bourdin possibly corrupting the file 8826823f3c5SBlaise Bourdin 883c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 8846823f3c5SBlaise Bourdin @*/ 885d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer) 886d71ae5a4SJacob Faibussowitsch { 8876823f3c5SBlaise Bourdin DM dm; 8886823f3c5SBlaise Bourdin MPI_Comm comm; 8896823f3c5SBlaise Bourdin PetscMPIInt rank; 8906823f3c5SBlaise Bourdin 8916823f3c5SBlaise Bourdin int exoid, offsetN = 0, offsetZ = 0; 8926823f3c5SBlaise Bourdin const char *vecname; 8936823f3c5SBlaise Bourdin PetscInt step; 8946823f3c5SBlaise Bourdin 8956823f3c5SBlaise Bourdin PetscFunctionBegin; 8969566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 8979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8989566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer, &exoid)); 8999566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 9009566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)v, &vecname)); 9016823f3c5SBlaise Bourdin 9029566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL)); 9039566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN)); 9049566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ)); 9051dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 9061dca8a05SBarry Smith if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN)); 9071dca8a05SBarry Smith else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ)); 9081dca8a05SBarry Smith else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 9093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 91078569bb4SMatthew G. Knepley } 91178569bb4SMatthew G. Knepley 91278569bb4SMatthew G. Knepley /* 91378569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 91478569bb4SMatthew G. Knepley 91520f4b53cSBarry Smith Collective 91678569bb4SMatthew G. Knepley 91778569bb4SMatthew G. Knepley Input Parameters: 91878569bb4SMatthew G. Knepley + v - The vector to be written 91978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9206823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 9216823f3c5SBlaise Bourdin - offset - the location of the variable in the file 92278569bb4SMatthew G. Knepley 92320f4b53cSBarry Smith Level: beginner 92420f4b53cSBarry Smith 92578569bb4SMatthew G. Knepley Notes: 92678569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 92778569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 92878569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 92978569bb4SMatthew G. Knepley amongst all nodal variables. 93078569bb4SMatthew G. Knepley 931c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 93278569bb4SMatthew G. Knepley @*/ 933d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 934d71ae5a4SJacob Faibussowitsch { 93578569bb4SMatthew G. Knepley MPI_Comm comm; 93678569bb4SMatthew G. Knepley PetscMPIInt size; 93778569bb4SMatthew G. Knepley DM dm; 93878569bb4SMatthew G. Knepley Vec vNatural, vComp; 93922a7544eSVaclav Hapla const PetscScalar *varray; 94078569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 94178569bb4SMatthew G. Knepley PetscBool useNatural; 94278569bb4SMatthew G. Knepley 94378569bb4SMatthew G. Knepley PetscFunctionBegin; 9449566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 9459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 9469566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 9479566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 94878569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 94978569bb4SMatthew G. Knepley if (useNatural) { 950b7352c5cSAlexis Marboeuf PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 9519566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 9529566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 95378569bb4SMatthew G. Knepley } else { 95478569bb4SMatthew G. Knepley vNatural = v; 95578569bb4SMatthew G. Knepley } 9566823f3c5SBlaise Bourdin 95778569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 95878569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 95978569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 9609566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 9619566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 96278569bb4SMatthew G. Knepley if (bs == 1) { 9639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 964792fecdfSBarry Smith PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 9659566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 96678569bb4SMatthew G. Knepley } else { 96778569bb4SMatthew G. Knepley IS compIS; 96878569bb4SMatthew G. Knepley PetscInt c; 96978569bb4SMatthew G. Knepley 9709566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 97178569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 9729566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 9739566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 9749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 975792fecdfSBarry Smith PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 9769566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 9779566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 97878569bb4SMatthew G. Knepley } 9799566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 98078569bb4SMatthew G. Knepley } 981b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(VecDestroy(&vNatural)); 9823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98378569bb4SMatthew G. Knepley } 98478569bb4SMatthew G. Knepley 98578569bb4SMatthew G. Knepley /* 98678569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 98778569bb4SMatthew G. Knepley 98820f4b53cSBarry Smith Collective 98978569bb4SMatthew G. Knepley 99078569bb4SMatthew G. Knepley Input Parameters: 99178569bb4SMatthew G. Knepley + v - The vector to be written 99278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9936823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 9946823f3c5SBlaise Bourdin - offset - the location of the variable in the file 99578569bb4SMatthew G. Knepley 99620f4b53cSBarry Smith Level: beginner 99720f4b53cSBarry Smith 99878569bb4SMatthew G. Knepley Notes: 99978569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 100078569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 100178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 100278569bb4SMatthew G. Knepley amongst all nodal variables. 100378569bb4SMatthew G. Knepley 1004db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 100578569bb4SMatthew G. Knepley */ 1006d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 1007d71ae5a4SJacob Faibussowitsch { 100878569bb4SMatthew G. Knepley MPI_Comm comm; 100978569bb4SMatthew G. Knepley PetscMPIInt size; 101078569bb4SMatthew G. Knepley DM dm; 101178569bb4SMatthew G. Knepley Vec vNatural, vComp; 101222a7544eSVaclav Hapla PetscScalar *varray; 101378569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 101478569bb4SMatthew G. Knepley PetscBool useNatural; 101578569bb4SMatthew G. Knepley 101678569bb4SMatthew G. Knepley PetscFunctionBegin; 10179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 10189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10199566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 10209566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 102178569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1022b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1023ad540459SPierre Jolivet else vNatural = v; 10246823f3c5SBlaise Bourdin 102578569bb4SMatthew G. Knepley /* Read local chunk from the file */ 10269566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 10279566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 102878569bb4SMatthew G. Knepley if (bs == 1) { 10299566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 1030792fecdfSBarry Smith PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray); 10319566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 103278569bb4SMatthew G. Knepley } else { 103378569bb4SMatthew G. Knepley IS compIS; 103478569bb4SMatthew G. Knepley PetscInt c; 103578569bb4SMatthew G. Knepley 10369566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 103778569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 10389566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 10399566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 10409566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 1041792fecdfSBarry Smith PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray); 10429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 10439566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 104478569bb4SMatthew G. Knepley } 10459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 104678569bb4SMatthew G. Knepley } 104778569bb4SMatthew G. Knepley if (useNatural) { 10489566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 10499566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1050b7352c5cSAlexis Marboeuf PetscCall(VecDestroy(&vNatural)); 105178569bb4SMatthew G. Knepley } 10523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105378569bb4SMatthew G. Knepley } 105478569bb4SMatthew G. Knepley 105578569bb4SMatthew G. Knepley /* 105678569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 105778569bb4SMatthew G. Knepley 105820f4b53cSBarry Smith Collective 105978569bb4SMatthew G. Knepley 106078569bb4SMatthew G. Knepley Input Parameters: 106178569bb4SMatthew G. Knepley + v - The vector to be written 106278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 10636823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 10646823f3c5SBlaise Bourdin - offset - the location of the variable in the file 106578569bb4SMatthew G. Knepley 106620f4b53cSBarry Smith Level: beginner 106720f4b53cSBarry Smith 106878569bb4SMatthew G. Knepley Notes: 106978569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 107078569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 107178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 107278569bb4SMatthew G. Knepley amongst all zonal variables. 107378569bb4SMatthew G. Knepley 1074c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 107578569bb4SMatthew G. Knepley */ 1076d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1077d71ae5a4SJacob Faibussowitsch { 107878569bb4SMatthew G. Knepley MPI_Comm comm; 107978569bb4SMatthew G. Knepley PetscMPIInt size; 108078569bb4SMatthew G. Knepley DM dm; 108178569bb4SMatthew G. Knepley Vec vNatural, vComp; 108222a7544eSVaclav Hapla const PetscScalar *varray; 108378569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 108478569bb4SMatthew G. Knepley PetscBool useNatural; 108578569bb4SMatthew G. Knepley IS compIS; 108678569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 108778569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 108878569bb4SMatthew G. Knepley 108978569bb4SMatthew G. Knepley PetscFunctionBegin; 10909566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 10919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10929566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 10939566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 109478569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 109578569bb4SMatthew G. Knepley if (useNatural) { 1096b7352c5cSAlexis Marboeuf PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 10979566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 10989566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 109978569bb4SMatthew G. Knepley } else { 110078569bb4SMatthew G. Knepley vNatural = v; 110178569bb4SMatthew G. Knepley } 11026823f3c5SBlaise Bourdin 110378569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 110478569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 110578569bb4SMatthew G. Knepley We assume that they are stored sequentially 110678569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1107a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 110878569bb4SMatthew G. Knepley to figure out what to save where. */ 110978569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 11109566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1111792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 111278569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 111378569bb4SMatthew G. Knepley ex_block block; 111478569bb4SMatthew G. Knepley 111578569bb4SMatthew G. Knepley block.id = csID[set]; 111678569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1117792fecdfSBarry Smith PetscCallExternal(ex_get_block_param, exoid, &block); 111878569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 111978569bb4SMatthew G. Knepley } 11209566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 11219566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 11229566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 112378569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 112478569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 112578569bb4SMatthew G. Knepley 112678569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 112778569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 112878569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 112978569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 113078569bb4SMatthew G. Knepley if (bs == 1) { 11319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 1132792fecdfSBarry 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)]); 11339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 113478569bb4SMatthew G. Knepley } else { 113578569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 11369566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 11379566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 11389566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 1139792fecdfSBarry 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)]); 11409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 11419566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 114278569bb4SMatthew G. Knepley } 114378569bb4SMatthew G. Knepley } 114478569bb4SMatthew G. Knepley csxs += csSize[set]; 114578569bb4SMatthew G. Knepley } 11469566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 11479566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 1148b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(VecDestroy(&vNatural)); 11493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115078569bb4SMatthew G. Knepley } 115178569bb4SMatthew G. Knepley 115278569bb4SMatthew G. Knepley /* 115378569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 115478569bb4SMatthew G. Knepley 115520f4b53cSBarry Smith Collective 115678569bb4SMatthew G. Knepley 115778569bb4SMatthew G. Knepley Input Parameters: 115878569bb4SMatthew G. Knepley + v - The vector to be written 115978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 11606823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 11616823f3c5SBlaise Bourdin - offset - the location of the variable in the file 116278569bb4SMatthew G. Knepley 116320f4b53cSBarry Smith Level: beginner 116420f4b53cSBarry Smith 116578569bb4SMatthew G. Knepley Notes: 116678569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 116778569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 116878569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 116978569bb4SMatthew G. Knepley amongst all zonal variables. 117078569bb4SMatthew G. Knepley 1171db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 117278569bb4SMatthew G. Knepley */ 1173d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 1174d71ae5a4SJacob Faibussowitsch { 117578569bb4SMatthew G. Knepley MPI_Comm comm; 117678569bb4SMatthew G. Knepley PetscMPIInt size; 117778569bb4SMatthew G. Knepley DM dm; 117878569bb4SMatthew G. Knepley Vec vNatural, vComp; 117922a7544eSVaclav Hapla PetscScalar *varray; 118078569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 118178569bb4SMatthew G. Knepley PetscBool useNatural; 118278569bb4SMatthew G. Knepley IS compIS; 118378569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 118478569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 118578569bb4SMatthew G. Knepley 118678569bb4SMatthew G. Knepley PetscFunctionBegin; 11879566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 11889566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11899566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 11909566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 119178569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 1192b7352c5cSAlexis Marboeuf if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural)); 1193ad540459SPierre Jolivet else vNatural = v; 11946823f3c5SBlaise Bourdin 119578569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 119678569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 119778569bb4SMatthew G. Knepley We assume that they are stored sequentially 119878569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1199a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 120078569bb4SMatthew G. Knepley to figure out what to save where. */ 120178569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 12029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1203792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID); 120478569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 120578569bb4SMatthew G. Knepley ex_block block; 120678569bb4SMatthew G. Knepley 120778569bb4SMatthew G. Knepley block.id = csID[set]; 120878569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1209792fecdfSBarry Smith PetscCallExternal(ex_get_block_param, exoid, &block); 121078569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 121178569bb4SMatthew G. Knepley } 12129566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 12139566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 12149566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS)); 121578569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 121678569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 121778569bb4SMatthew G. Knepley 121878569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 121978569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 122078569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 122178569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs)); 122278569bb4SMatthew G. Knepley if (bs == 1) { 12239566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 1224792fecdfSBarry 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)]); 12259566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 122678569bb4SMatthew G. Knepley } else { 122778569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 12289566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs)); 12299566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 12309566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 1231792fecdfSBarry 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)]); 12329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 12339566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 123478569bb4SMatthew G. Knepley } 123578569bb4SMatthew G. Knepley } 123678569bb4SMatthew G. Knepley csxs += csSize[set]; 123778569bb4SMatthew G. Knepley } 12389566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 12399566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 124078569bb4SMatthew G. Knepley if (useNatural) { 12419566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 12429566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 1243b7352c5cSAlexis Marboeuf PetscCall(VecDestroy(&vNatural)); 124478569bb4SMatthew G. Knepley } 12453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 124678569bb4SMatthew G. Knepley } 1247b53c8456SSatish Balay #endif 124878569bb4SMatthew G. Knepley 12491e50132fSMatthew G. Knepley /*@ 1250a1cb98faSBarry Smith PetscViewerExodusIIGetId - Get the file id of the `PETSCVIEWEREXODUSII` file 12511e50132fSMatthew G. Knepley 125220f4b53cSBarry Smith Logically Collective 12531e50132fSMatthew G. Knepley 12541e50132fSMatthew G. Knepley Input Parameter: 1255a1cb98faSBarry Smith . viewer - the `PetscViewer` 12561e50132fSMatthew G. Knepley 12571e50132fSMatthew G. Knepley Output Parameter: 1258d8d19677SJose E. Roman . exoid - The ExodusII file id 12591e50132fSMatthew G. Knepley 12601e50132fSMatthew G. Knepley Level: intermediate 12611e50132fSMatthew G. Knepley 1262a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 12631e50132fSMatthew G. Knepley @*/ 1264d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) 1265d71ae5a4SJacob Faibussowitsch { 12661e50132fSMatthew G. Knepley PetscFunctionBegin; 12671e50132fSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1268cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetId_C", (PetscViewer, int *), (viewer, exoid)); 12693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12706823f3c5SBlaise Bourdin } 12716823f3c5SBlaise Bourdin 12726823f3c5SBlaise Bourdin /*@ 12736823f3c5SBlaise Bourdin PetscViewerExodusIISetOrder - Set the elements order in the exodusII file. 12746823f3c5SBlaise Bourdin 12756823f3c5SBlaise Bourdin Collective 12766823f3c5SBlaise Bourdin 12776823f3c5SBlaise Bourdin Input Parameters: 1278a1cb98faSBarry Smith + viewer - the `PETSCVIEWEREXODUSII` viewer 127998a6a78aSPierre Jolivet - order - elements order 12806823f3c5SBlaise Bourdin 12816823f3c5SBlaise Bourdin Output Parameter: 12826823f3c5SBlaise Bourdin 12836823f3c5SBlaise Bourdin Level: beginner 12846823f3c5SBlaise Bourdin 1285*42747ad1SJacob Faibussowitsch .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()` 12866823f3c5SBlaise Bourdin @*/ 1287d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order) 1288d71ae5a4SJacob Faibussowitsch { 12896823f3c5SBlaise Bourdin PetscFunctionBegin; 12906823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1291cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerSetOrder_C", (PetscViewer, PetscInt), (viewer, order)); 12923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12936823f3c5SBlaise Bourdin } 12946823f3c5SBlaise Bourdin 12956823f3c5SBlaise Bourdin /*@ 12966823f3c5SBlaise Bourdin PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file. 12976823f3c5SBlaise Bourdin 12986823f3c5SBlaise Bourdin Collective 12996823f3c5SBlaise Bourdin 13006823f3c5SBlaise Bourdin Input Parameters: 1301a1cb98faSBarry Smith + viewer - the `PETSCVIEWEREXODUSII` viewer 130298a6a78aSPierre Jolivet - order - elements order 13036823f3c5SBlaise Bourdin 13046823f3c5SBlaise Bourdin Output Parameter: 13056823f3c5SBlaise Bourdin 13066823f3c5SBlaise Bourdin Level: beginner 13076823f3c5SBlaise Bourdin 1308*42747ad1SJacob Faibussowitsch .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIISetOrder()` 13096823f3c5SBlaise Bourdin @*/ 1310d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order) 1311d71ae5a4SJacob Faibussowitsch { 13126823f3c5SBlaise Bourdin PetscFunctionBegin; 13136823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1314cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetOrder_C", (PetscViewer, PetscInt *), (viewer, order)); 13153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13166823f3c5SBlaise Bourdin } 13176823f3c5SBlaise Bourdin 13186823f3c5SBlaise Bourdin /*@C 13196823f3c5SBlaise Bourdin PetscViewerExodusIIOpen - Opens a file for ExodusII input/output. 13206823f3c5SBlaise Bourdin 13216823f3c5SBlaise Bourdin Collective 13226823f3c5SBlaise Bourdin 13236823f3c5SBlaise Bourdin Input Parameters: 13246823f3c5SBlaise Bourdin + comm - MPI communicator 13256823f3c5SBlaise Bourdin . name - name of file 13266823f3c5SBlaise Bourdin - type - type of file 1327a1cb98faSBarry Smith .vb 1328a1cb98faSBarry Smith FILE_MODE_WRITE - create new file for binary output 1329a1cb98faSBarry Smith FILE_MODE_READ - open existing file for binary input 1330a1cb98faSBarry Smith FILE_MODE_APPEND - open existing file for binary output 1331a1cb98faSBarry Smith .ve 13326823f3c5SBlaise Bourdin 13336823f3c5SBlaise Bourdin Output Parameter: 1334a1cb98faSBarry Smith . exo - `PETSCVIEWEREXODUSII` `PetscViewer` for Exodus II input/output to use with the specified file 13356823f3c5SBlaise Bourdin 13366823f3c5SBlaise Bourdin Level: beginner 13376823f3c5SBlaise Bourdin 1338a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 133960225df5SJacob Faibussowitsch `DMLoad()`, `PetscFileMode`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 13406823f3c5SBlaise Bourdin @*/ 1341d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) 1342d71ae5a4SJacob Faibussowitsch { 13436823f3c5SBlaise Bourdin PetscFunctionBegin; 13449566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, exo)); 13459566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII)); 13469566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*exo, type)); 13479566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*exo, name)); 13489566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*exo)); 13493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13501e50132fSMatthew G. Knepley } 13511e50132fSMatthew G. Knepley 135233751fbdSMichael Lange /*@C 1353a1cb98faSBarry Smith DMPlexCreateExodusFromFile - Create a `DMPLEX` mesh from an ExodusII file. 135433751fbdSMichael Lange 1355d083f849SBarry Smith Collective 135633751fbdSMichael Lange 135733751fbdSMichael Lange Input Parameters: 135833751fbdSMichael Lange + comm - The MPI communicator 135933751fbdSMichael Lange . filename - The name of the ExodusII file 136033751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 136133751fbdSMichael Lange 136233751fbdSMichael Lange Output Parameter: 1363a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 136433751fbdSMichael Lange 136533751fbdSMichael Lange Level: beginner 136633751fbdSMichael Lange 13671cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()` 136833751fbdSMichael Lange @*/ 1369d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 1370d71ae5a4SJacob Faibussowitsch { 137133751fbdSMichael Lange PetscMPIInt rank; 137233751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1373e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 137433751fbdSMichael Lange float version; 137533751fbdSMichael Lange #endif 137633751fbdSMichael Lange 137733751fbdSMichael Lange PetscFunctionBegin; 13784f572ea9SToby Isaac PetscAssertPointer(filename, 2); 13799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 138033751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1381dd400576SPatrick Sanan if (rank == 0) { 138233751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 138308401ef6SPierre Jolivet PetscCheck(exoid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 138433751fbdSMichael Lange } 13859566063dSJacob Faibussowitsch PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm)); 138648a46eb9SPierre Jolivet if (rank == 0) PetscCallExternal(ex_close, exoid); 13873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 138833751fbdSMichael Lange #else 138933751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 139033751fbdSMichael Lange #endif 139133751fbdSMichael Lange } 139233751fbdSMichael Lange 13938f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 1394d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct) 1395d71ae5a4SJacob Faibussowitsch { 13968f861fbcSMatthew G. Knepley PetscBool flg; 13978f861fbcSMatthew G. Knepley 13988f861fbcSMatthew G. Knepley PetscFunctionBegin; 13998f861fbcSMatthew G. Knepley *ct = DM_POLYTOPE_UNKNOWN; 1400207ab81fSDavid Andrs PetscCall(PetscStrcmp(elem_type, "BAR2", &flg)); 1401207ab81fSDavid Andrs if (flg) { 1402207ab81fSDavid Andrs *ct = DM_POLYTOPE_SEGMENT; 1403207ab81fSDavid Andrs goto done; 1404207ab81fSDavid Andrs } 1405207ab81fSDavid Andrs PetscCall(PetscStrcmp(elem_type, "BAR3", &flg)); 1406207ab81fSDavid Andrs if (flg) { 1407207ab81fSDavid Andrs *ct = DM_POLYTOPE_SEGMENT; 1408207ab81fSDavid Andrs goto done; 1409207ab81fSDavid Andrs } 14109566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI", &flg)); 14119371c9d4SSatish Balay if (flg) { 14129371c9d4SSatish Balay *ct = DM_POLYTOPE_TRIANGLE; 14139371c9d4SSatish Balay goto done; 14149371c9d4SSatish Balay } 14159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI3", &flg)); 14169371c9d4SSatish Balay if (flg) { 14179371c9d4SSatish Balay *ct = DM_POLYTOPE_TRIANGLE; 14189371c9d4SSatish Balay goto done; 14199371c9d4SSatish Balay } 14209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD", &flg)); 14219371c9d4SSatish Balay if (flg) { 14229371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 14239371c9d4SSatish Balay goto done; 14249371c9d4SSatish Balay } 14259566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg)); 14269371c9d4SSatish Balay if (flg) { 14279371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 14289371c9d4SSatish Balay goto done; 14299371c9d4SSatish Balay } 14309566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg)); 14319371c9d4SSatish Balay if (flg) { 14329371c9d4SSatish Balay *ct = DM_POLYTOPE_QUADRILATERAL; 14339371c9d4SSatish Balay goto done; 14349371c9d4SSatish Balay } 14359566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TETRA", &flg)); 14369371c9d4SSatish Balay if (flg) { 14379371c9d4SSatish Balay *ct = DM_POLYTOPE_TETRAHEDRON; 14389371c9d4SSatish Balay goto done; 14399371c9d4SSatish Balay } 14409566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TET4", &flg)); 14419371c9d4SSatish Balay if (flg) { 14429371c9d4SSatish Balay *ct = DM_POLYTOPE_TETRAHEDRON; 14439371c9d4SSatish Balay goto done; 14449371c9d4SSatish Balay } 14459566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg)); 14469371c9d4SSatish Balay if (flg) { 14479371c9d4SSatish Balay *ct = DM_POLYTOPE_TRI_PRISM; 14489371c9d4SSatish Balay goto done; 14499371c9d4SSatish Balay } 14509566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX", &flg)); 14519371c9d4SSatish Balay if (flg) { 14529371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14539371c9d4SSatish Balay goto done; 14549371c9d4SSatish Balay } 14559566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX8", &flg)); 14569371c9d4SSatish Balay if (flg) { 14579371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14589371c9d4SSatish Balay goto done; 14599371c9d4SSatish Balay } 14609566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg)); 14619371c9d4SSatish Balay if (flg) { 14629371c9d4SSatish Balay *ct = DM_POLYTOPE_HEXAHEDRON; 14639371c9d4SSatish Balay goto done; 14649371c9d4SSatish Balay } 146598921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type); 14668f861fbcSMatthew G. Knepley done: 14673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14688f861fbcSMatthew G. Knepley } 14698f861fbcSMatthew G. Knepley #endif 14708f861fbcSMatthew G. Knepley 1471552f7358SJed Brown /*@ 1472a1cb98faSBarry Smith DMPlexCreateExodus - Create a `DMPLEX` mesh from an ExodusII file ID. 1473552f7358SJed Brown 1474d083f849SBarry Smith Collective 1475552f7358SJed Brown 1476552f7358SJed Brown Input Parameters: 1477552f7358SJed Brown + comm - The MPI communicator 1478552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 1479552f7358SJed Brown - interpolate - Create faces and edges in the mesh 1480552f7358SJed Brown 1481552f7358SJed Brown Output Parameter: 1482a1cb98faSBarry Smith . dm - The `DM` object representing the mesh 1483552f7358SJed Brown 1484552f7358SJed Brown Level: beginner 1485552f7358SJed Brown 148660225df5SJacob Faibussowitsch .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMCreate()` 1487552f7358SJed Brown @*/ 1488d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 1489d71ae5a4SJacob Faibussowitsch { 1490552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 1491552f7358SJed Brown PetscMPIInt num_proc, rank; 1492ae9eebabSLisandro Dalcin DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL; 1493552f7358SJed Brown PetscSection coordSection; 1494552f7358SJed Brown Vec coordinates; 1495552f7358SJed Brown PetscScalar *coords; 1496552f7358SJed Brown PetscInt coordSize, v; 1497552f7358SJed Brown /* Read from ex_get_init() */ 1498552f7358SJed Brown char title[PETSC_MAX_PATH_LEN + 1]; 14995f80ce2aSJacob Faibussowitsch int dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0; 1500552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 1501552f7358SJed Brown #endif 1502552f7358SJed Brown 1503552f7358SJed Brown PetscFunctionBegin; 1504552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 15059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 15069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 15079566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 15089566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1509a5b23f4aSJose E. Roman /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */ 1510dd400576SPatrick Sanan if (rank == 0) { 15119566063dSJacob Faibussowitsch PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1)); 1512792fecdfSBarry Smith PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 151328b400f6SJacob Faibussowitsch PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set"); 1514552f7358SJed Brown } 15159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm)); 15169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 15179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)*dm, title)); 15189566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices)); 1519412e9a14SMatthew G. Knepley /* We do not want this label automatically computed, instead we compute it here */ 15209566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "celltype")); 1521552f7358SJed Brown 1522552f7358SJed Brown /* Read cell sets information */ 1523dd400576SPatrick Sanan if (rank == 0) { 1524552f7358SJed Brown PetscInt *cone; 1525e8f6893fSMatthew G. Knepley int c, cs, ncs, c_loc, v, v_loc; 1526552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 1527e8f6893fSMatthew G. Knepley int *cs_id, *cs_order; 1528552f7358SJed Brown /* Read from ex_get_elem_block() */ 1529552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN + 1]; 1530e8f6893fSMatthew G. Knepley int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr; 1531552f7358SJed Brown /* Read from ex_get_elem_conn() */ 1532552f7358SJed Brown int *cs_connect; 1533552f7358SJed Brown 1534552f7358SJed Brown /* Get cell sets IDs */ 15359566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order)); 1536792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id); 1537552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 1538552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 1539e8f6893fSMatthew G. Knepley /* Check for a hybrid mesh */ 1540e8f6893fSMatthew G. Knepley for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) { 15418f861fbcSMatthew G. Knepley DMPolytopeType ct; 15428f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 15438f861fbcSMatthew G. Knepley 15449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1545792fecdfSBarry Smith PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 15469566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 15478f861fbcSMatthew G. Knepley dim = PetscMax(dim, DMPolytopeTypeGetDim(ct)); 1548792fecdfSBarry 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); 15498f861fbcSMatthew G. Knepley switch (ct) { 15508f861fbcSMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 1551e8f6893fSMatthew G. Knepley cs_order[cs] = cs; 1552e8f6893fSMatthew G. Knepley ++num_hybrid; 1553e8f6893fSMatthew G. Knepley break; 1554e8f6893fSMatthew G. Knepley default: 1555e8f6893fSMatthew G. Knepley for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1]; 1556e8f6893fSMatthew G. Knepley cs_order[cs - num_hybrid] = cs; 1557e8f6893fSMatthew G. Knepley } 1558e8f6893fSMatthew G. Knepley } 1559552f7358SJed Brown /* First set sizes */ 1560e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 15618f861fbcSMatthew G. Knepley DMPolytopeType ct; 15628f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 1563e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 15648f861fbcSMatthew G. Knepley 15659566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1566792fecdfSBarry Smith PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type); 15679566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 1568792fecdfSBarry 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); 1569552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 15709566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell)); 15719566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, c, ct)); 1572552f7358SJed Brown } 1573552f7358SJed Brown } 15749566063dSJacob Faibussowitsch for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 15759566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 1576e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1577e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 1578792fecdfSBarry 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); 15799566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone)); 1580792fecdfSBarry Smith PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL); 1581eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 1582552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1583636e64ffSStefano Zampini DMPolytopeType ct; 1584636e64ffSStefano Zampini 1585ad540459SPierre Jolivet for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1; 15869566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(*dm, c, &ct)); 15879566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCell(ct, cone)); 15889566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, cone)); 15899566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs])); 1590552f7358SJed Brown } 15919566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_connect, cone)); 1592552f7358SJed Brown } 15939566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_id, cs_order)); 1594552f7358SJed Brown } 15958f861fbcSMatthew G. Knepley { 15968f861fbcSMatthew G. Knepley PetscInt ints[] = {dim, dimEmbed}; 15978f861fbcSMatthew G. Knepley 15989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm)); 15999566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, ints[0])); 16009566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, ints[1])); 16018f861fbcSMatthew G. Knepley dim = ints[0]; 16028f861fbcSMatthew G. Knepley dimEmbed = ints[1]; 16038f861fbcSMatthew G. Knepley } 16049566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 16059566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1606552f7358SJed Brown if (interpolate) { 16075fd9971aSMatthew G. Knepley DM idm; 1608552f7358SJed Brown 16099566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 16109566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1611552f7358SJed Brown *dm = idm; 1612552f7358SJed Brown } 1613552f7358SJed Brown 1614552f7358SJed Brown /* Create vertex set label */ 1615dd400576SPatrick Sanan if (rank == 0 && (num_vs > 0)) { 1616552f7358SJed Brown int vs, v; 1617552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 1618552f7358SJed Brown int *vs_id; 1619552f7358SJed Brown /* Read from ex_get_node_set_param() */ 1620f958083aSBlaise Bourdin int num_vertex_in_set; 1621552f7358SJed Brown /* Read from ex_get_node_set() */ 1622552f7358SJed Brown int *vs_vertex_list; 1623552f7358SJed Brown 1624552f7358SJed Brown /* Get vertex set ids */ 16259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vs, &vs_id)); 1626792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id); 1627552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 1628792fecdfSBarry Smith PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 16299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list)); 1630792fecdfSBarry Smith PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL); 163148a46eb9SPierre 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])); 16329566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_vertex_list)); 1633552f7358SJed Brown } 16349566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_id)); 1635552f7358SJed Brown } 1636552f7358SJed Brown /* Read coordinates */ 16379566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 16389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 16399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed)); 16409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 1641552f7358SJed Brown for (v = numCells; v < numCells + numVertices; ++v) { 16429566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed)); 16439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed)); 1644552f7358SJed Brown } 16459566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 16469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 16479566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 16489566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates")); 16499566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 16509566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, dimEmbed)); 16519566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates, VECSTANDARD)); 16529566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 1653dd400576SPatrick Sanan if (rank == 0) { 1654e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 1655552f7358SJed Brown 16569566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z)); 1657792fecdfSBarry Smith PetscCallExternal(ex_get_coord, exoid, x, y, z); 16588f861fbcSMatthew G. Knepley if (dimEmbed > 0) { 16598f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v]; 16600d644c17SKarl Rupp } 16618f861fbcSMatthew G. Knepley if (dimEmbed > 1) { 16628f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v]; 16630d644c17SKarl Rupp } 16648f861fbcSMatthew G. Knepley if (dimEmbed > 2) { 16658f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v]; 16660d644c17SKarl Rupp } 16679566063dSJacob Faibussowitsch PetscCall(PetscFree3(x, y, z)); 1668552f7358SJed Brown } 16699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 16709566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 16719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1672552f7358SJed Brown 1673552f7358SJed Brown /* Create side set label */ 1674dd400576SPatrick Sanan if (rank == 0 && interpolate && (num_fs > 0)) { 1675552f7358SJed Brown int fs, f, voff; 1676552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 1677552f7358SJed Brown int *fs_id; 1678552f7358SJed Brown /* Read from ex_get_side_set_param() */ 1679f958083aSBlaise Bourdin int num_side_in_set; 1680552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 1681552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 1682ef073283Smichael_afanasiev /* Read side set labels */ 1683ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH + 1]; 1684034d25ccSMatthew G. Knepley size_t fs_name_len; 1685552f7358SJed Brown 1686552f7358SJed Brown /* Get side set ids */ 16879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_fs, &fs_id)); 1688792fecdfSBarry Smith PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id); 168982cad2ffSMatthew G. Knepley // Ids 1 and 2 are reserved by ExodusII for indicating things in 3D 1690552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 1691792fecdfSBarry Smith PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL); 16929566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list)); 1693792fecdfSBarry Smith PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list); 1694ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 1695ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1696034d25ccSMatthew G. Knepley if (!fs_name_err) { 1697034d25ccSMatthew G. Knepley PetscCall(PetscStrlen(fs_name, &fs_name_len)); 1698034d25ccSMatthew G. Knepley if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH)); 1699034d25ccSMatthew G. Knepley } 170082cad2ffSMatthew G. Knepley PetscCheck(fs_id[fs] != 1 && fs_id[fs] != 2, comm, PETSC_ERR_ARG_WRONG, "Side set %s marker cannot be %d since this is reserved by ExodusII", fs_name, fs_id[fs]); 1701552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 17020298fd71SBarry Smith const PetscInt *faces = NULL; 1703552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1704552f7358SJed Brown PetscInt faceVertices[4], v; 1705552f7358SJed Brown 1706197f6770SJed Brown PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize); 1707ad540459SPierre Jolivet for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1; 17089566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1709197f6770SJed Brown PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces); 17109566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs])); 1711ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 1712034d25ccSMatthew G. Knepley if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs])); 17139566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1714552f7358SJed Brown } 17159566063dSJacob Faibussowitsch PetscCall(PetscFree2(fs_vertex_count_list, fs_vertex_list)); 1716552f7358SJed Brown } 17179566063dSJacob Faibussowitsch PetscCall(PetscFree(fs_id)); 1718552f7358SJed Brown } 1719ae9eebabSLisandro Dalcin 1720ae9eebabSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 17219371c9d4SSatish Balay enum { 17229371c9d4SSatish Balay n = 3 17239371c9d4SSatish Balay }; 1724ae9eebabSLisandro Dalcin PetscBool flag[n]; 1725ae9eebabSLisandro Dalcin 1726ae9eebabSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 1727ae9eebabSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 1728ae9eebabSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 17299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 17309566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 17319566063dSJacob Faibussowitsch if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 17329566063dSJacob Faibussowitsch if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 1733ae9eebabSLisandro Dalcin } 17343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1735552f7358SJed Brown #else 1736552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1737552f7358SJed Brown #endif 1738552f7358SJed Brown } 1739