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) 1278569bb4SMatthew G. Knepley /* 131e50132fSMatthew G. Knepley PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator. 141e50132fSMatthew G. Knepley 151e50132fSMatthew G. Knepley Collective 161e50132fSMatthew G. Knepley 171e50132fSMatthew G. Knepley Input Parameter: 181e50132fSMatthew G. Knepley . comm - the MPI communicator to share the ExodusII PetscViewer 191e50132fSMatthew G. Knepley 201e50132fSMatthew G. Knepley Level: intermediate 211e50132fSMatthew G. Knepley 221e50132fSMatthew G. Knepley Notes: 231e50132fSMatthew G. Knepley misses Fortran bindings 241e50132fSMatthew G. Knepley 251e50132fSMatthew G. Knepley Notes: 261e50132fSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return 271e50132fSMatthew G. Knepley an error code. The GLVIS PetscViewer is usually used in the form 281e50132fSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm)); 291e50132fSMatthew G. Knepley 30db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()` 311e50132fSMatthew G. Knepley */ 321e50132fSMatthew G. Knepley PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm) 331e50132fSMatthew G. Knepley { 341e50132fSMatthew G. Knepley PetscViewer viewer; 351e50132fSMatthew G. Knepley PetscErrorCode ierr; 361e50132fSMatthew G. Knepley 371e50132fSMatthew G. Knepley PetscFunctionBegin; 381e50132fSMatthew G. Knepley ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer); 3900373969SVaclav Hapla if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 401e50132fSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject) viewer); 4100373969SVaclav Hapla if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 421e50132fSMatthew G. Knepley PetscFunctionReturn(viewer); 431e50132fSMatthew G. Knepley } 441e50132fSMatthew G. Knepley 451e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer) 461e50132fSMatthew G. Knepley { 476823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) v->data; 481e50132fSMatthew G. Knepley 491e50132fSMatthew G. Knepley PetscFunctionBegin; 509566063dSJacob Faibussowitsch if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename)); 519566063dSJacob Faibussowitsch if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid: %d\n", exo->exoid)); 529566063dSJacob Faibussowitsch if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode: %d\n", exo->btype)); 539566063dSJacob Faibussowitsch if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order: %d\n", exo->order)); 541e50132fSMatthew G. Knepley PetscFunctionReturn(0); 551e50132fSMatthew G. Knepley } 561e50132fSMatthew G. Knepley 571e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v) 581e50132fSMatthew G. Knepley { 591e50132fSMatthew G. Knepley PetscFunctionBegin; 60d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options"); 61d0609cedSBarry Smith PetscOptionsHeadEnd(); 621e50132fSMatthew G. Knepley PetscFunctionReturn(0); 631e50132fSMatthew G. Knepley } 641e50132fSMatthew G. Knepley 651e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer) 661e50132fSMatthew G. Knepley { 671e50132fSMatthew G. Knepley PetscFunctionBegin; 681e50132fSMatthew G. Knepley PetscFunctionReturn(0); 691e50132fSMatthew G. Knepley } 701e50132fSMatthew G. Knepley 711e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer) 721e50132fSMatthew G. Knepley { 731e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 741e50132fSMatthew G. Knepley 751e50132fSMatthew G. Knepley PetscFunctionBegin; 76a74df02fSJacob Faibussowitsch if (exo->exoid >= 0) {PetscStackCallStandard(ex_close,exo->exoid);} 779566063dSJacob Faibussowitsch PetscCall(PetscFree(exo->filename)); 789566063dSJacob Faibussowitsch PetscCall(PetscFree(exo)); 799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL)); 809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL)); 819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL)); 82*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL)); 839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetId_C",NULL)); 849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetOrder_C",NULL)); 859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerSetOrder_C",NULL)); 861e50132fSMatthew G. Knepley PetscFunctionReturn(0); 871e50132fSMatthew G. Knepley } 881e50132fSMatthew G. Knepley 891e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[]) 901e50132fSMatthew G. Knepley { 911e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 921e50132fSMatthew G. Knepley PetscMPIInt rank; 931e50132fSMatthew G. Knepley int CPU_word_size, IO_word_size, EXO_mode; 946823f3c5SBlaise Bourdin MPI_Info mpi_info = MPI_INFO_NULL; 956823f3c5SBlaise Bourdin float EXO_version; 961e50132fSMatthew G. Knepley 979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank)); 981e50132fSMatthew G. Knepley CPU_word_size = sizeof(PetscReal); 991e50132fSMatthew G. Knepley IO_word_size = sizeof(PetscReal); 1001e50132fSMatthew G. Knepley 1011e50132fSMatthew G. Knepley PetscFunctionBegin; 1026823f3c5SBlaise Bourdin if (exo->exoid >= 0) { 103a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_close,exo->exoid); 1046823f3c5SBlaise Bourdin exo->exoid = -1; 1056823f3c5SBlaise Bourdin } 1069566063dSJacob Faibussowitsch if (exo->filename) PetscCall(PetscFree(exo->filename)); 1079566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &exo->filename)); 1081e50132fSMatthew G. Knepley switch (exo->btype) { 1091e50132fSMatthew G. Knepley case FILE_MODE_READ: 1106823f3c5SBlaise Bourdin EXO_mode = EX_READ; 1111e50132fSMatthew G. Knepley break; 1121e50132fSMatthew G. Knepley case FILE_MODE_APPEND: 1136823f3c5SBlaise Bourdin case FILE_MODE_UPDATE: 1146823f3c5SBlaise Bourdin case FILE_MODE_APPEND_UPDATE: 1156823f3c5SBlaise Bourdin /* Will fail if the file does not already exist */ 1166823f3c5SBlaise Bourdin EXO_mode = EX_WRITE; 1171e50132fSMatthew G. Knepley break; 1181e50132fSMatthew G. Knepley case FILE_MODE_WRITE: 1196823f3c5SBlaise Bourdin /* 1206823f3c5SBlaise Bourdin exodus only allows writing geometry upon file creation, so we will let DMView create the file. 1216823f3c5SBlaise Bourdin */ 1226823f3c5SBlaise Bourdin PetscFunctionReturn(0); 1231e50132fSMatthew G. Knepley break; 1241e50132fSMatthew G. Knepley default: 1251e50132fSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 1261e50132fSMatthew G. Knepley } 1271e50132fSMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 1281e50132fSMatthew G. Knepley EXO_mode += EX_ALL_INT64_API; 1291e50132fSMatthew G. Knepley #endif 1306823f3c5SBlaise Bourdin exo->exoid = ex_open_par(name,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,PETSC_COMM_WORLD,mpi_info); 13108401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name); 1321e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1331e50132fSMatthew G. Knepley } 1341e50132fSMatthew G. Knepley 1351e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name) 1361e50132fSMatthew G. Knepley { 1371e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1381e50132fSMatthew G. Knepley 1391e50132fSMatthew G. Knepley PetscFunctionBegin; 1401e50132fSMatthew G. Knepley *name = exo->filename; 1411e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1421e50132fSMatthew G. Knepley } 1431e50132fSMatthew G. Knepley 1441e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type) 1451e50132fSMatthew G. Knepley { 1461e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1471e50132fSMatthew G. Knepley 1481e50132fSMatthew G. Knepley PetscFunctionBegin; 1491e50132fSMatthew G. Knepley exo->btype = type; 1501e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1511e50132fSMatthew G. Knepley } 1521e50132fSMatthew G. Knepley 1531e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type) 1541e50132fSMatthew G. Knepley { 1551e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1561e50132fSMatthew G. Knepley 1571e50132fSMatthew G. Knepley PetscFunctionBegin; 1581e50132fSMatthew G. Knepley *type = exo->btype; 1591e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1601e50132fSMatthew G. Knepley } 1611e50132fSMatthew G. Knepley 1621e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid) 1631e50132fSMatthew G. Knepley { 1641e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1651e50132fSMatthew G. Knepley 1661e50132fSMatthew G. Knepley PetscFunctionBegin; 1671e50132fSMatthew G. Knepley *exoid = exo->exoid; 1681e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1691e50132fSMatthew G. Knepley } 1701e50132fSMatthew G. Knepley 1716823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order) 1721e50132fSMatthew G. Knepley { 1736823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1741e50132fSMatthew G. Knepley 1751e50132fSMatthew G. Knepley PetscFunctionBegin; 1766823f3c5SBlaise Bourdin *order = exo->order; 1776823f3c5SBlaise Bourdin PetscFunctionReturn(0); 1786823f3c5SBlaise Bourdin } 1796823f3c5SBlaise Bourdin 1806823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order) 1816823f3c5SBlaise Bourdin { 1826823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1836823f3c5SBlaise Bourdin 1846823f3c5SBlaise Bourdin PetscFunctionBegin; 1856823f3c5SBlaise Bourdin exo->order = order; 1861e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1871e50132fSMatthew G. Knepley } 1881e50132fSMatthew G. Knepley 1891e50132fSMatthew G. Knepley /*MC 19000373969SVaclav Hapla PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file 1911e50132fSMatthew G. Knepley 192db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`, 193db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 1941e50132fSMatthew G. Knepley 1951e50132fSMatthew G. Knepley Level: beginner 1961e50132fSMatthew G. Knepley M*/ 1971e50132fSMatthew G. Knepley 1981e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v) 1991e50132fSMatthew G. Knepley { 2001e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo; 2011e50132fSMatthew G. Knepley 2021e50132fSMatthew G. Knepley PetscFunctionBegin; 2039566063dSJacob Faibussowitsch PetscCall(PetscNewLog(v,&exo)); 2041e50132fSMatthew G. Knepley 2051e50132fSMatthew G. Knepley v->data = (void*) exo; 2061e50132fSMatthew G. Knepley v->ops->destroy = PetscViewerDestroy_ExodusII; 2071e50132fSMatthew G. Knepley v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII; 2081e50132fSMatthew G. Knepley v->ops->setup = PetscViewerSetUp_ExodusII; 2091e50132fSMatthew G. Knepley v->ops->view = PetscViewerView_ExodusII; 2101e50132fSMatthew G. Knepley v->ops->flush = 0; 2111e50132fSMatthew G. Knepley exo->btype = (PetscFileMode) -1; 2121e50132fSMatthew G. Knepley exo->filename = 0; 2131e50132fSMatthew G. Knepley exo->exoid = -1; 2141e50132fSMatthew G. Knepley 2159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII)); 2169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII)); 2179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII)); 2189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII)); 2199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetId_C",PetscViewerExodusIIGetId_ExodusII)); 2209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerSetOrder_C",PetscViewerExodusIISetOrder_ExodusII)); 2219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetOrder_C",PetscViewerExodusIIGetOrder_ExodusII)); 2221e50132fSMatthew G. Knepley PetscFunctionReturn(0); 2231e50132fSMatthew G. Knepley } 2241e50132fSMatthew G. Knepley 2251e50132fSMatthew G. Knepley /* 22678569bb4SMatthew G. Knepley EXOGetVarIndex - Locate a result in an exodus file based on its name 22778569bb4SMatthew G. Knepley 2286823f3c5SBlaise Bourdin Collective 22978569bb4SMatthew G. Knepley 23078569bb4SMatthew G. Knepley Input Parameters: 23178569bb4SMatthew G. Knepley + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 23278569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 23378569bb4SMatthew G. Knepley - name - the name of the result 23478569bb4SMatthew G. Knepley 23578569bb4SMatthew G. Knepley Output Parameters: 23678569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result 23778569bb4SMatthew G. Knepley 23878569bb4SMatthew G. Knepley Notes: 23978569bb4SMatthew G. Knepley The exodus variable index is obtained by comparing name and the 24078569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance if name is "V" 24178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 24278569bb4SMatthew G. Knepley amongst all variables of type obj_type. 24378569bb4SMatthew G. Knepley 24478569bb4SMatthew G. Knepley Level: beginner 24578569bb4SMatthew G. Knepley 246c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 24778569bb4SMatthew G. Knepley */ 2486823f3c5SBlaise Bourdin PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 24978569bb4SMatthew G. Knepley { 2506de834b4SMatthew G. Knepley int num_vars, i, j; 25178569bb4SMatthew G. Knepley char ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1]; 25278569bb4SMatthew G. Knepley const int num_suffix = 5; 25378569bb4SMatthew G. Knepley char *suffix[5]; 25478569bb4SMatthew G. Knepley PetscBool flg; 25578569bb4SMatthew G. Knepley 25678569bb4SMatthew G. Knepley PetscFunctionBegin; 25778569bb4SMatthew G. Knepley suffix[0] = (char *) ""; 25878569bb4SMatthew G. Knepley suffix[1] = (char *) "_X"; 25978569bb4SMatthew G. Knepley suffix[2] = (char *) "_XX"; 26078569bb4SMatthew G. Knepley suffix[3] = (char *) "_1"; 26178569bb4SMatthew G. Knepley suffix[4] = (char *) "_11"; 26278569bb4SMatthew G. Knepley 2636823f3c5SBlaise Bourdin *varIndex = -1; 264a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_variable_param,exoid, obj_type, &num_vars); 26578569bb4SMatthew G. Knepley for (i = 0; i < num_vars; ++i) { 266a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_variable_name,exoid, obj_type, i+1, var_name); 26778569bb4SMatthew G. Knepley for (j = 0; j < num_suffix; ++j) { 2689566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH)); 2699566063dSJacob Faibussowitsch PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH)); 2709566063dSJacob Faibussowitsch PetscCall(PetscStrcasecmp(ext_name, var_name, &flg)); 27178569bb4SMatthew G. Knepley if (flg) { 27278569bb4SMatthew G. Knepley *varIndex = i+1; 2736823f3c5SBlaise Bourdin } 2746823f3c5SBlaise Bourdin } 2756823f3c5SBlaise Bourdin } 27678569bb4SMatthew G. Knepley PetscFunctionReturn(0); 27778569bb4SMatthew G. Knepley } 27878569bb4SMatthew G. Knepley 27978569bb4SMatthew G. Knepley /* 2806823f3c5SBlaise Bourdin DMView_PlexExodusII - Write a DM to disk in exodus format 28178569bb4SMatthew G. Knepley 28278569bb4SMatthew G. Knepley Collective on dm 28378569bb4SMatthew G. Knepley 28478569bb4SMatthew G. Knepley Input Parameters: 28578569bb4SMatthew G. Knepley + dm - The dm to be written 2866823f3c5SBlaise Bourdin . viewer - an exodusII viewer 28778569bb4SMatthew G. Knepley 28878569bb4SMatthew G. Knepley Notes: 28978569bb4SMatthew G. Knepley Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 29078569bb4SMatthew G. Knepley consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 29178569bb4SMatthew G. Knepley 2926823f3c5SBlaise Bourdin If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices) 2936823f3c5SBlaise Bourdin will be written. 29478569bb4SMatthew G. Knepley 29578569bb4SMatthew G. Knepley DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 29678569bb4SMatthew G. Knepley on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 2976823f3c5SBlaise Bourdin The order of the mesh shall be set using PetscViewerExodusIISetOrder 29878569bb4SMatthew G. Knepley It should be extended to use PetscFE objects. 29978569bb4SMatthew G. Knepley 3006823f3c5SBlaise Bourdin This function will only handle TRI, TET, QUAD, and HEX cells. 30178569bb4SMatthew G. Knepley Level: beginner 30278569bb4SMatthew G. Knepley 3036823f3c5SBlaise Bourdin .seealso: 30478569bb4SMatthew G. Knepley */ 3056823f3c5SBlaise Bourdin PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer) 30678569bb4SMatthew G. Knepley { 30778569bb4SMatthew G. Knepley enum ElemType {TRI, QUAD, TET, HEX}; 30878569bb4SMatthew G. Knepley MPI_Comm comm; 3096823f3c5SBlaise Bourdin PetscInt degree; /* the order of the mesh */ 31078569bb4SMatthew G. Knepley /* Connectivity Variables */ 31178569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 31278569bb4SMatthew G. Knepley /* Cell Sets */ 31378569bb4SMatthew G. Knepley DMLabel csLabel; 31478569bb4SMatthew G. Knepley IS csIS; 31578569bb4SMatthew G. Knepley const PetscInt *csIdx; 31678569bb4SMatthew G. Knepley PetscInt num_cs, cs; 31778569bb4SMatthew G. Knepley enum ElemType *type; 318fe209ef9SBlaise Bourdin PetscBool hasLabel; 31978569bb4SMatthew G. Knepley /* Coordinate Variables */ 32078569bb4SMatthew G. Knepley DM cdm; 3216823f3c5SBlaise Bourdin PetscSection coordSection; 32278569bb4SMatthew G. Knepley Vec coord; 32378569bb4SMatthew G. Knepley PetscInt **nodes; 324eae8a387SMatthew G. Knepley PetscInt depth, d, dim, skipCells = 0; 32578569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 32678569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 32778569bb4SMatthew G. Knepley PetscMPIInt rank, size; 32878569bb4SMatthew G. Knepley const char *dmName; 329fe209ef9SBlaise Bourdin PetscInt nodesTriP1[4] = {3,0,0,0}; 330fe209ef9SBlaise Bourdin PetscInt nodesTriP2[4] = {3,3,0,0}; 331fe209ef9SBlaise Bourdin PetscInt nodesQuadP1[4] = {4,0,0,0}; 332fe209ef9SBlaise Bourdin PetscInt nodesQuadP2[4] = {4,4,0,1}; 333fe209ef9SBlaise Bourdin PetscInt nodesTetP1[4] = {4,0,0,0}; 334fe209ef9SBlaise Bourdin PetscInt nodesTetP2[4] = {4,6,0,0}; 335fe209ef9SBlaise Bourdin PetscInt nodesHexP1[4] = {8,0,0,0}; 336fe209ef9SBlaise Bourdin PetscInt nodesHexP2[4] = {8,12,6,1}; 3376823f3c5SBlaise Bourdin int CPU_word_size, IO_word_size, EXO_mode; 3386823f3c5SBlaise Bourdin float EXO_version; 3396823f3c5SBlaise Bourdin 3406823f3c5SBlaise Bourdin PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 3416823f3c5SBlaise Bourdin 34278569bb4SMatthew G. Knepley PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 3449566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 3466823f3c5SBlaise Bourdin 3476823f3c5SBlaise Bourdin /* 3486823f3c5SBlaise Bourdin Creating coordSection is a collective operation so we do it somewhat out of sequence 3496823f3c5SBlaise Bourdin */ 3509566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &coordSection)); 3519566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalSetUp(dm)); 352c5853193SPierre Jolivet if (rank == 0) { 3536823f3c5SBlaise Bourdin switch (exo->btype) { 3546823f3c5SBlaise Bourdin case FILE_MODE_READ: 3556823f3c5SBlaise Bourdin case FILE_MODE_APPEND: 3566823f3c5SBlaise Bourdin case FILE_MODE_UPDATE: 3576823f3c5SBlaise Bourdin case FILE_MODE_APPEND_UPDATE: 3586823f3c5SBlaise Bourdin /* exodusII does not allow writing geometry to an existing file */ 35998921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename); 3606823f3c5SBlaise Bourdin case FILE_MODE_WRITE: 3616823f3c5SBlaise Bourdin /* Create an empty file if one already exists*/ 3626823f3c5SBlaise Bourdin EXO_mode = EX_CLOBBER; 3636823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 3646823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 3656823f3c5SBlaise Bourdin #endif 3666823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 3676823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 3686823f3c5SBlaise Bourdin exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size); 36908401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename); 3706823f3c5SBlaise Bourdin 3716823f3c5SBlaise Bourdin break; 3726823f3c5SBlaise Bourdin default: 3736823f3c5SBlaise Bourdin SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3746823f3c5SBlaise Bourdin } 3756823f3c5SBlaise Bourdin 37678569bb4SMatthew G. Knepley /* --- Get DM info --- */ 3779566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) dm, &dmName)); 3789566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 3799566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 3809566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 3819566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 3829566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd)); 3839566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd)); 3849566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd)); 38578569bb4SMatthew G. Knepley numCells = cEnd - cStart; 38678569bb4SMatthew G. Knepley numEdges = eEnd - eStart; 38778569bb4SMatthew G. Knepley numVertices = vEnd - vStart; 38878569bb4SMatthew G. Knepley if (depth == 3) {numFaces = fEnd - fStart;} 38978569bb4SMatthew G. Knepley else {numFaces = 0;} 3909566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs)); 3919566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs)); 3929566063dSJacob Faibussowitsch PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs)); 3939566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coord)); 3949566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 39578569bb4SMatthew G. Knepley if (num_cs > 0) { 3969566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel)); 3979566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(csLabel, &csIS)); 3989566063dSJacob Faibussowitsch PetscCall(ISGetIndices(csIS, &csIdx)); 39978569bb4SMatthew G. Knepley } 4009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &nodes)); 40178569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 4029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_cs, &type)); 40378569bb4SMatthew G. Knepley numNodes = numVertices; 4046823f3c5SBlaise Bourdin 4059566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetOrder(viewer, °ree)); 40678569bb4SMatthew G. Knepley if (degree == 2) {numNodes += numEdges;} 40778569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 40878569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 40978569bb4SMatthew G. Knepley IS stratumIS; 41078569bb4SMatthew G. Knepley const PetscInt *cells; 41178569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 41278569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 41378569bb4SMatthew G. Knepley 4149566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 4159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 4169566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 4179566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 41878569bb4SMatthew G. Knepley switch (dim) { 41978569bb4SMatthew G. Knepley case 2: 42078569bb4SMatthew G. Knepley if (closureSize == 3*dim) {type[cs] = TRI;} 42178569bb4SMatthew G. Knepley else if (closureSize == 4*dim) {type[cs] = QUAD;} 42263a3b9bcSJacob Faibussowitsch 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); 42378569bb4SMatthew G. Knepley break; 42478569bb4SMatthew G. Knepley case 3: 42578569bb4SMatthew G. Knepley if (closureSize == 4*dim) {type[cs] = TET;} 42678569bb4SMatthew G. Knepley else if (closureSize == 8*dim) {type[cs] = HEX;} 42763a3b9bcSJacob Faibussowitsch 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); 42878569bb4SMatthew G. Knepley break; 42963a3b9bcSJacob Faibussowitsch default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim); 43078569bb4SMatthew G. Knepley } 43178569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;} 43278569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == HEX)) {numNodes += csSize; numNodes += numFaces;} 4339566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz)); 43478569bb4SMatthew G. Knepley /* Set nodes and Element type */ 43578569bb4SMatthew G. Knepley if (type[cs] == TRI) { 43678569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 43778569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 43878569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 43978569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 44078569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 44178569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 44278569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 44378569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 44478569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 44578569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 44678569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 44778569bb4SMatthew G. Knepley } 44878569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 44978569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3]*csSize; 45078569bb4SMatthew G. Knepley 4519566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 4529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 45378569bb4SMatthew G. Knepley } 454a74df02fSJacob Faibussowitsch if (num_cs > 0) {PetscStackCallStandard(ex_put_init,exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);} 45578569bb4SMatthew G. Knepley /* --- Connectivity --- */ 45678569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 45778569bb4SMatthew G. Knepley IS stratumIS; 45878569bb4SMatthew G. Knepley const PetscInt *cells; 459fe209ef9SBlaise Bourdin PetscInt *connect, off = 0; 460eae8a387SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 46178569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 46278569bb4SMatthew G. Knepley char *elem_type = NULL; 46378569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 46478569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 46578569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 46678569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 46778569bb4SMatthew G. Knepley 4689566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 4699566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 4709566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 47178569bb4SMatthew G. Knepley /* Set Element type */ 47278569bb4SMatthew G. Knepley if (type[cs] == TRI) { 47378569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 47478569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 47578569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 47678569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 47778569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 47878569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 47978569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 48078569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 48178569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 48278569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 48378569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 48478569bb4SMatthew G. Knepley } 48578569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 4869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect)); 487a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_block,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1); 48878569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 48978569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 49078569bb4SMatthew G. Knepley if (depth > 1) { 49178569bb4SMatthew G. Knepley if (dim == 2) { 4929566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure)); 49378569bb4SMatthew G. Knepley } else if (dim == 3) { 49478569bb4SMatthew G. Knepley PetscInt *closure = NULL; 49578569bb4SMatthew G. Knepley 4969566063dSJacob Faibussowitsch PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure)); 4979566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 49878569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 4999566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure)); 50078569bb4SMatthew G. Knepley } 50178569bb4SMatthew G. Knepley } 50278569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 50378569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 50478569bb4SMatthew G. Knepley PetscInt *closure = NULL; 50578569bb4SMatthew G. Knepley PetscInt temp, i; 50678569bb4SMatthew G. Knepley 5079566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 50878569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 50978569bb4SMatthew G. Knepley if (i < nodes[cs][0]) {/* Vertices */ 510fe209ef9SBlaise Bourdin connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 511fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 51278569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 513fe209ef9SBlaise Bourdin connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 514fe209ef9SBlaise Bourdin if (nodes[cs][2] == 0) connect[i+off] -= numFaces; 515fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 51678569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]) { /* Cells */ 517fe209ef9SBlaise Bourdin connect[i+off] = closure[0] + 1; 518fe209ef9SBlaise Bourdin connect[i+off] -= skipCells; 51978569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]) { /* Faces */ 520fe209ef9SBlaise Bourdin connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 521fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 52278569bb4SMatthew G. Knepley } else { 523fe209ef9SBlaise Bourdin connect[i+off] = -1; 52478569bb4SMatthew G. Knepley } 52578569bb4SMatthew G. Knepley } 52678569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 52778569bb4SMatthew G. Knepley if (type[cs] == TET) { 528fe209ef9SBlaise Bourdin temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp; 52978569bb4SMatthew G. Knepley if (degree == 2) { 530e2c9586dSBlaise Bourdin temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp; 531e2c9586dSBlaise Bourdin temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp; 53278569bb4SMatthew G. Knepley } 53378569bb4SMatthew G. Knepley } 53478569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 53578569bb4SMatthew G. Knepley if (type[cs] == HEX) { 536fe209ef9SBlaise Bourdin temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp; 53778569bb4SMatthew G. Knepley if (degree == 2) { 538fe209ef9SBlaise Bourdin temp = connect[8+off]; connect[8+off] = connect[11+off]; connect[11+off] = temp; 539fe209ef9SBlaise Bourdin temp = connect[9+off]; connect[9+off] = connect[10+off]; connect[10+off] = temp; 540fe209ef9SBlaise Bourdin temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp; 541fe209ef9SBlaise Bourdin temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp; 54278569bb4SMatthew G. Knepley 543fe209ef9SBlaise Bourdin temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp; 544fe209ef9SBlaise Bourdin temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp; 545fe209ef9SBlaise Bourdin temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp; 546fe209ef9SBlaise Bourdin temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp; 54778569bb4SMatthew G. Knepley 548fe209ef9SBlaise Bourdin temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp; 549fe209ef9SBlaise Bourdin temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp; 550fe209ef9SBlaise Bourdin temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp; 55178569bb4SMatthew G. Knepley } 55278569bb4SMatthew G. Knepley } 553fe209ef9SBlaise Bourdin off += connectSize; 5549566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure)); 55578569bb4SMatthew G. Knepley } 556a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_conn,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0); 55778569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0)*csSize; 5589566063dSJacob Faibussowitsch PetscCall(PetscFree(connect)); 5599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 5609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 56178569bb4SMatthew G. Knepley } 5629566063dSJacob Faibussowitsch PetscCall(PetscFree(type)); 56378569bb4SMatthew G. Knepley /* --- Coordinates --- */ 5649566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd)); 5651e50132fSMatthew G. Knepley if (num_cs) { 56678569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 5679566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd)); 56878569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 5699566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0)); 57078569bb4SMatthew G. Knepley } 57178569bb4SMatthew G. Knepley } 5721e50132fSMatthew G. Knepley } 57378569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 57478569bb4SMatthew G. Knepley IS stratumIS; 57578569bb4SMatthew G. Knepley const PetscInt *cells; 57678569bb4SMatthew G. Knepley PetscInt csSize, c; 57778569bb4SMatthew G. Knepley 5789566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS)); 5799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &cells)); 5809566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &csSize)); 58178569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 5829566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0)); 58378569bb4SMatthew G. Knepley } 5849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &cells)); 5859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 58678569bb4SMatthew G. Knepley } 58778569bb4SMatthew G. Knepley if (num_cs > 0) { 5889566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(csIS, &csIdx)); 5899566063dSJacob Faibussowitsch PetscCall(ISDestroy(&csIS)); 59078569bb4SMatthew G. Knepley } 5919566063dSJacob Faibussowitsch PetscCall(PetscFree(nodes)); 5929566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 59378569bb4SMatthew G. Knepley if (numNodes > 0) { 59478569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 595233c95e0SFrancesco Ballarin PetscScalar *closure, *cval; 596233c95e0SFrancesco Ballarin PetscReal *coords; 59778569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 59878569bb4SMatthew G. Knepley 5996823f3c5SBlaise Bourdin /* There can't be more than 24 values in the closure of a point for the coord coordSection */ 6009566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure)); 6019566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord)); 6029566063dSJacob Faibussowitsch PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 60378569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 6049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(coordSection, p, &hasDof)); 60578569bb4SMatthew G. Knepley if (hasDof) { 60678569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 60778569bb4SMatthew G. Knepley 6089566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure)); 60978569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 61078569bb4SMatthew G. Knepley cval[d] = 0.0; 61178569bb4SMatthew G. Knepley for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 612233c95e0SFrancesco Ballarin coords[d*numNodes+n] = PetscRealPart(cval[d]) * dim / closureSize; 61378569bb4SMatthew G. Knepley } 61478569bb4SMatthew G. Knepley ++n; 61578569bb4SMatthew G. Knepley } 61678569bb4SMatthew G. Knepley } 617a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_coord,exo->exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]); 6189566063dSJacob Faibussowitsch PetscCall(PetscFree3(coords, cval, closure)); 619a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_coord_names,exo->exoid, (char **) coordNames); 62078569bb4SMatthew G. Knepley } 6216823f3c5SBlaise Bourdin 622fe209ef9SBlaise Bourdin /* --- Node Sets/Vertex Sets --- */ 623fe209ef9SBlaise Bourdin DMHasLabel(dm, "Vertex Sets", &hasLabel); 624fe209ef9SBlaise Bourdin if (hasLabel) { 625fe209ef9SBlaise Bourdin PetscInt i, vs, vsSize; 626fe209ef9SBlaise Bourdin const PetscInt *vsIdx, *vertices; 627fe209ef9SBlaise Bourdin PetscInt *nodeList; 628fe209ef9SBlaise Bourdin IS vsIS, stratumIS; 629fe209ef9SBlaise Bourdin DMLabel vsLabel; 6309566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel)); 6319566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(vsLabel, &vsIS)); 6329566063dSJacob Faibussowitsch PetscCall(ISGetIndices(vsIS, &vsIdx)); 633fe209ef9SBlaise Bourdin for (vs=0; vs<num_vs; ++vs) { 6349566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS)); 6359566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &vertices)); 6369566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &vsSize)); 6379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(vsSize, &nodeList)); 638fe209ef9SBlaise Bourdin for (i=0; i<vsSize; ++i) { 639fe209ef9SBlaise Bourdin nodeList[i] = vertices[i] - skipCells + 1; 640fe209ef9SBlaise Bourdin } 641a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0); 642a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_set,exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL); 6439566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &vertices)); 6449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 6459566063dSJacob Faibussowitsch PetscCall(PetscFree(nodeList)); 646fe209ef9SBlaise Bourdin } 6479566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(vsIS, &vsIdx)); 6489566063dSJacob Faibussowitsch PetscCall(ISDestroy(&vsIS)); 649fe209ef9SBlaise Bourdin } 650fe209ef9SBlaise Bourdin /* --- Side Sets/Face Sets --- */ 6519566063dSJacob Faibussowitsch PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel)); 652fe209ef9SBlaise Bourdin if (hasLabel) { 653fe209ef9SBlaise Bourdin PetscInt i, j, fs, fsSize; 654fe209ef9SBlaise Bourdin const PetscInt *fsIdx, *faces; 655fe209ef9SBlaise Bourdin IS fsIS, stratumIS; 656fe209ef9SBlaise Bourdin DMLabel fsLabel; 657fe209ef9SBlaise Bourdin PetscInt numPoints, *points; 658fe209ef9SBlaise Bourdin PetscInt elem_list_size = 0; 659fe209ef9SBlaise Bourdin PetscInt *elem_list, *elem_ind, *side_list; 660fe209ef9SBlaise Bourdin 6619566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel)); 662fe209ef9SBlaise Bourdin /* Compute size of Node List and Element List */ 6639566063dSJacob Faibussowitsch PetscCall(DMLabelGetValueIS(fsLabel, &fsIS)); 6649566063dSJacob Faibussowitsch PetscCall(ISGetIndices(fsIS, &fsIdx)); 665fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 6669566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 6679566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 668fe209ef9SBlaise Bourdin elem_list_size += fsSize; 6699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 670fe209ef9SBlaise Bourdin } 671d0609cedSBarry Smith PetscCall(PetscMalloc3(num_fs, &elem_ind,elem_list_size, &elem_list,elem_list_size, &side_list)); 672fe209ef9SBlaise Bourdin elem_ind[0] = 0; 673fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 6749566063dSJacob Faibussowitsch PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS)); 6759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(stratumIS, &faces)); 6769566063dSJacob Faibussowitsch PetscCall(ISGetSize(stratumIS, &fsSize)); 677fe209ef9SBlaise Bourdin /* Set Parameters */ 678a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0); 679fe209ef9SBlaise Bourdin /* Indices */ 680fe209ef9SBlaise Bourdin if (fs<num_fs-1) { 681fe209ef9SBlaise Bourdin elem_ind[fs+1] = elem_ind[fs] + fsSize; 682fe209ef9SBlaise Bourdin } 683fe209ef9SBlaise Bourdin 684fe209ef9SBlaise Bourdin for (i=0; i<fsSize; ++i) { 685fe209ef9SBlaise Bourdin /* Element List */ 686fe209ef9SBlaise Bourdin points = NULL; 6879566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 688fe209ef9SBlaise Bourdin elem_list[elem_ind[fs] + i] = points[2] +1; 6899566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points)); 690fe209ef9SBlaise Bourdin 691fe209ef9SBlaise Bourdin /* Side List */ 692fe209ef9SBlaise Bourdin points = NULL; 6939566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points)); 694fe209ef9SBlaise Bourdin for (j=1; j<numPoints; ++j) { 695fe209ef9SBlaise Bourdin if (points[j*2]==faces[i]) {break;} 696fe209ef9SBlaise Bourdin } 697fe209ef9SBlaise Bourdin /* Convert HEX sides */ 698fe209ef9SBlaise Bourdin if (numPoints == 27) { 699fe209ef9SBlaise Bourdin if (j == 1) {j = 5;} 700fe209ef9SBlaise Bourdin else if (j == 2) {j = 6;} 701fe209ef9SBlaise Bourdin else if (j == 3) {j = 1;} 702fe209ef9SBlaise Bourdin else if (j == 4) {j = 3;} 703fe209ef9SBlaise Bourdin else if (j == 5) {j = 2;} 704fe209ef9SBlaise Bourdin else if (j == 6) {j = 4;} 705fe209ef9SBlaise Bourdin } 706fe209ef9SBlaise Bourdin /* Convert TET sides */ 707fe209ef9SBlaise Bourdin if (numPoints == 15) { 708fe209ef9SBlaise Bourdin --j; 709fe209ef9SBlaise Bourdin if (j == 0) {j = 4;} 710fe209ef9SBlaise Bourdin } 711fe209ef9SBlaise Bourdin side_list[elem_ind[fs] + i] = j; 7129566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points)); 713fe209ef9SBlaise Bourdin 714fe209ef9SBlaise Bourdin } 7159566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(stratumIS, &faces)); 7169566063dSJacob Faibussowitsch PetscCall(ISDestroy(&stratumIS)); 717fe209ef9SBlaise Bourdin } 7189566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(fsIS, &fsIdx)); 7199566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fsIS)); 720fe209ef9SBlaise Bourdin 721fe209ef9SBlaise Bourdin /* Put side sets */ 722fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 723a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_set,exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]); 724fe209ef9SBlaise Bourdin } 7259566063dSJacob Faibussowitsch PetscCall(PetscFree3(elem_ind,elem_list,side_list)); 726fe209ef9SBlaise Bourdin } 7276823f3c5SBlaise Bourdin /* 7286823f3c5SBlaise Bourdin close the exodus file 7296823f3c5SBlaise Bourdin */ 7306823f3c5SBlaise Bourdin ex_close(exo->exoid); 7316823f3c5SBlaise Bourdin exo->exoid = -1; 7326823f3c5SBlaise Bourdin } 7339566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&coordSection)); 7346823f3c5SBlaise Bourdin 7356823f3c5SBlaise Bourdin /* 7366823f3c5SBlaise Bourdin reopen the file in parallel 7376823f3c5SBlaise Bourdin */ 7386823f3c5SBlaise Bourdin EXO_mode = EX_WRITE; 7396823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES) 7406823f3c5SBlaise Bourdin EXO_mode += EX_ALL_INT64_API; 7416823f3c5SBlaise Bourdin #endif 7426823f3c5SBlaise Bourdin CPU_word_size = sizeof(PetscReal); 7436823f3c5SBlaise Bourdin IO_word_size = sizeof(PetscReal); 7446823f3c5SBlaise Bourdin exo->exoid = ex_open_par(exo->filename,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,comm,MPI_INFO_NULL); 74508401ef6SPierre Jolivet PetscCheck(exo->exoid >= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename); 7466823f3c5SBlaise Bourdin PetscFunctionReturn(0); 7476823f3c5SBlaise Bourdin } 7486823f3c5SBlaise Bourdin 7496823f3c5SBlaise Bourdin /* 7506823f3c5SBlaise Bourdin VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 7516823f3c5SBlaise Bourdin 7526823f3c5SBlaise Bourdin Collective on v 7536823f3c5SBlaise Bourdin 7546823f3c5SBlaise Bourdin Input Parameters: 7556823f3c5SBlaise Bourdin + v - The vector to be written 7566823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 7576823f3c5SBlaise Bourdin 7586823f3c5SBlaise Bourdin Notes: 7596823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 7606823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 7616823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 7626823f3c5SBlaise Bourdin amongst all variables. 7636823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 7646823f3c5SBlaise Bourdin possibly corrupting the file 7656823f3c5SBlaise Bourdin 7666823f3c5SBlaise Bourdin Level: beginner 7676823f3c5SBlaise Bourdin 768c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 7696823f3c5SBlaise Bourdin @*/ 7706823f3c5SBlaise Bourdin PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer) 7716823f3c5SBlaise Bourdin { 7726823f3c5SBlaise Bourdin DM dm; 7736823f3c5SBlaise Bourdin MPI_Comm comm; 7746823f3c5SBlaise Bourdin PetscMPIInt rank; 7756823f3c5SBlaise Bourdin 7766823f3c5SBlaise Bourdin int exoid,offsetN = 0, offsetZ = 0; 7776823f3c5SBlaise Bourdin const char *vecname; 7786823f3c5SBlaise Bourdin PetscInt step; 7796823f3c5SBlaise Bourdin 7806823f3c5SBlaise Bourdin PetscFunctionBegin; 7819566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) v, &comm)); 7829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 7839566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer,&exoid)); 7849566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 7859566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) v, &vecname)); 7866823f3c5SBlaise Bourdin 7879566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL)); 7889566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN)); 7899566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ)); 7901dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 7916823f3c5SBlaise Bourdin if (offsetN > 0) { 7929566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN)); 7936823f3c5SBlaise Bourdin } else if (offsetZ > 0) { 7949566063dSJacob Faibussowitsch PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ)); 79598921bdaSJacob Faibussowitsch } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 7966823f3c5SBlaise Bourdin PetscFunctionReturn(0); 7976823f3c5SBlaise Bourdin } 7986823f3c5SBlaise Bourdin 7996823f3c5SBlaise Bourdin /* 8006823f3c5SBlaise Bourdin VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file 8016823f3c5SBlaise Bourdin 8026823f3c5SBlaise Bourdin Collective on v 8036823f3c5SBlaise Bourdin 8046823f3c5SBlaise Bourdin Input Parameters: 8056823f3c5SBlaise Bourdin + v - The vector to be written 8066823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file 8076823f3c5SBlaise Bourdin 8086823f3c5SBlaise Bourdin Notes: 8096823f3c5SBlaise Bourdin The exodus result variable index is obtained by comparing the Vec name and the 8106823f3c5SBlaise Bourdin names of variables declared in the exodus file. For instance for a Vec named "V" 8116823f3c5SBlaise Bourdin the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 8126823f3c5SBlaise Bourdin amongst all variables. 8136823f3c5SBlaise Bourdin In the event where a nodal and zonal variable both match, the function will return an error instead of 8146823f3c5SBlaise Bourdin possibly corrupting the file 8156823f3c5SBlaise Bourdin 8166823f3c5SBlaise Bourdin Level: beginner 8176823f3c5SBlaise Bourdin 818c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()` 8196823f3c5SBlaise Bourdin @*/ 8206823f3c5SBlaise Bourdin PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer) 8216823f3c5SBlaise Bourdin { 8226823f3c5SBlaise Bourdin DM dm; 8236823f3c5SBlaise Bourdin MPI_Comm comm; 8246823f3c5SBlaise Bourdin PetscMPIInt rank; 8256823f3c5SBlaise Bourdin 8266823f3c5SBlaise Bourdin int exoid,offsetN = 0, offsetZ = 0; 8276823f3c5SBlaise Bourdin const char *vecname; 8286823f3c5SBlaise Bourdin PetscInt step; 8296823f3c5SBlaise Bourdin 8306823f3c5SBlaise Bourdin PetscFunctionBegin; 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) v, &comm)); 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 8339566063dSJacob Faibussowitsch PetscCall(PetscViewerExodusIIGetId(viewer,&exoid)); 8349566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8359566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject) v, &vecname)); 8366823f3c5SBlaise Bourdin 8379566063dSJacob Faibussowitsch PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL)); 8389566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN)); 8399566063dSJacob Faibussowitsch PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ)); 8401dca8a05SBarry Smith PetscCheck(offsetN > 0 || offsetZ > 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname); 8411dca8a05SBarry Smith if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN)); 8421dca8a05SBarry Smith else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ)); 8431dca8a05SBarry Smith else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname); 84478569bb4SMatthew G. Knepley PetscFunctionReturn(0); 84578569bb4SMatthew G. Knepley } 84678569bb4SMatthew G. Knepley 84778569bb4SMatthew G. Knepley /* 84878569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 84978569bb4SMatthew G. Knepley 85078569bb4SMatthew G. Knepley Collective on v 85178569bb4SMatthew G. Knepley 85278569bb4SMatthew G. Knepley Input Parameters: 85378569bb4SMatthew G. Knepley + v - The vector to be written 85478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 8556823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 8566823f3c5SBlaise Bourdin - offset - the location of the variable in the file 85778569bb4SMatthew G. Knepley 85878569bb4SMatthew G. Knepley Notes: 85978569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 86078569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 86178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 86278569bb4SMatthew G. Knepley amongst all nodal variables. 86378569bb4SMatthew G. Knepley 86478569bb4SMatthew G. Knepley Level: beginner 86578569bb4SMatthew G. Knepley 866c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 86778569bb4SMatthew G. Knepley @*/ 8686823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 86978569bb4SMatthew G. Knepley { 87078569bb4SMatthew G. Knepley MPI_Comm comm; 87178569bb4SMatthew G. Knepley PetscMPIInt size; 87278569bb4SMatthew G. Knepley DM dm; 87378569bb4SMatthew G. Knepley Vec vNatural, vComp; 87422a7544eSVaclav Hapla const PetscScalar *varray; 87578569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 87678569bb4SMatthew G. Knepley PetscBool useNatural; 87778569bb4SMatthew G. Knepley 87878569bb4SMatthew G. Knepley PetscFunctionBegin; 8799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) v, &comm)); 8809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 8819566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 8829566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 88378569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 88478569bb4SMatthew G. Knepley if (useNatural) { 8859566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &vNatural)); 8869566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 8879566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 88878569bb4SMatthew G. Knepley } else { 88978569bb4SMatthew G. Knepley vNatural = v; 89078569bb4SMatthew G. Knepley } 8916823f3c5SBlaise Bourdin 89278569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 89378569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 89478569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 8959566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 8969566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 89778569bb4SMatthew G. Knepley if (bs == 1) { 8989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 899a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray); 9009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 90178569bb4SMatthew G. Knepley } else { 90278569bb4SMatthew G. Knepley IS compIS; 90378569bb4SMatthew G. Knepley PetscInt c; 90478569bb4SMatthew G. Knepley 9059566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS)); 90678569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 9079566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs)); 9089566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 9099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 910a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray); 9119566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 9129566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 91378569bb4SMatthew G. Knepley } 9149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 91578569bb4SMatthew G. Knepley } 9169566063dSJacob Faibussowitsch if (useNatural) PetscCall(DMRestoreGlobalVector(dm, &vNatural)); 91778569bb4SMatthew G. Knepley PetscFunctionReturn(0); 91878569bb4SMatthew G. Knepley } 91978569bb4SMatthew G. Knepley 92078569bb4SMatthew G. Knepley /* 92178569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 92278569bb4SMatthew G. Knepley 92378569bb4SMatthew G. Knepley Collective on v 92478569bb4SMatthew G. Knepley 92578569bb4SMatthew G. Knepley Input Parameters: 92678569bb4SMatthew G. Knepley + v - The vector to be written 92778569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9286823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 9296823f3c5SBlaise Bourdin - offset - the location of the variable in the file 93078569bb4SMatthew G. Knepley 93178569bb4SMatthew G. Knepley Notes: 93278569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 93378569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 93478569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 93578569bb4SMatthew G. Knepley amongst all nodal variables. 93678569bb4SMatthew G. Knepley 93778569bb4SMatthew G. Knepley Level: beginner 93878569bb4SMatthew G. Knepley 939db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()` 94078569bb4SMatthew G. Knepley */ 9416823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) 94278569bb4SMatthew G. Knepley { 94378569bb4SMatthew G. Knepley MPI_Comm comm; 94478569bb4SMatthew G. Knepley PetscMPIInt size; 94578569bb4SMatthew G. Knepley DM dm; 94678569bb4SMatthew G. Knepley Vec vNatural, vComp; 94722a7544eSVaclav Hapla PetscScalar *varray; 94878569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 94978569bb4SMatthew G. Knepley PetscBool useNatural; 95078569bb4SMatthew G. Knepley 95178569bb4SMatthew G. Knepley PetscFunctionBegin; 9529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) v, &comm)); 9539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 9549566063dSJacob Faibussowitsch PetscCall(VecGetDM(v,&dm)); 9559566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 95678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 9579566063dSJacob Faibussowitsch if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural)); 95878569bb4SMatthew G. Knepley else {vNatural = v;} 9596823f3c5SBlaise Bourdin 96078569bb4SMatthew G. Knepley /* Read local chunk from the file */ 9619566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 9629566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 96378569bb4SMatthew G. Knepley if (bs == 1) { 9649566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 965a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray); 9669566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 96778569bb4SMatthew G. Knepley } else { 96878569bb4SMatthew G. Knepley IS compIS; 96978569bb4SMatthew G. Knepley PetscInt c; 97078569bb4SMatthew G. Knepley 9719566063dSJacob Faibussowitsch PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS)); 97278569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 9739566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs)); 9749566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 9759566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 976a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray); 9779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 9789566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 97978569bb4SMatthew G. Knepley } 9809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&compIS)); 98178569bb4SMatthew G. Knepley } 98278569bb4SMatthew G. Knepley if (useNatural) { 9839566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 9849566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 9859566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &vNatural)); 98678569bb4SMatthew G. Knepley } 98778569bb4SMatthew G. Knepley PetscFunctionReturn(0); 98878569bb4SMatthew G. Knepley } 98978569bb4SMatthew G. Knepley 99078569bb4SMatthew G. Knepley /* 99178569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 99278569bb4SMatthew G. Knepley 99378569bb4SMatthew G. Knepley Collective on v 99478569bb4SMatthew G. Knepley 99578569bb4SMatthew G. Knepley Input Parameters: 99678569bb4SMatthew G. Knepley + v - The vector to be written 99778569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 9986823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1) 9996823f3c5SBlaise Bourdin - offset - the location of the variable in the file 100078569bb4SMatthew G. Knepley 100178569bb4SMatthew G. Knepley Notes: 100278569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 100378569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 100478569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 100578569bb4SMatthew G. Knepley amongst all zonal variables. 100678569bb4SMatthew G. Knepley 100778569bb4SMatthew G. Knepley Level: beginner 100878569bb4SMatthew G. Knepley 1009c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 101078569bb4SMatthew G. Knepley */ 10116823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 101278569bb4SMatthew G. Knepley { 101378569bb4SMatthew G. Knepley MPI_Comm comm; 101478569bb4SMatthew G. Knepley PetscMPIInt size; 101578569bb4SMatthew G. Knepley DM dm; 101678569bb4SMatthew G. Knepley Vec vNatural, vComp; 101722a7544eSVaclav Hapla const PetscScalar *varray; 101878569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 101978569bb4SMatthew G. Knepley PetscBool useNatural; 102078569bb4SMatthew G. Knepley IS compIS; 102178569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 102278569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 102378569bb4SMatthew G. Knepley 102478569bb4SMatthew G. Knepley PetscFunctionBegin; 10259566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v, &comm)); 10269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 10279566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 10289566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 102978569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 103078569bb4SMatthew G. Knepley if (useNatural) { 10319566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &vNatural)); 10329566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural)); 10339566063dSJacob Faibussowitsch PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural)); 103478569bb4SMatthew G. Knepley } else { 103578569bb4SMatthew G. Knepley vNatural = v; 103678569bb4SMatthew G. Knepley } 10376823f3c5SBlaise Bourdin 103878569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 103978569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 104078569bb4SMatthew G. Knepley We assume that they are stored sequentially 104178569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1042a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 104378569bb4SMatthew G. Knepley to figure out what to save where. */ 104478569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 10459566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1046a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID); 104778569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 104878569bb4SMatthew G. Knepley ex_block block; 104978569bb4SMatthew G. Knepley 105078569bb4SMatthew G. Knepley block.id = csID[set]; 105178569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1052a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_block_param,exoid, &block); 105378569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 105478569bb4SMatthew G. Knepley } 10559566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 10569566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 10579566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS)); 105878569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 105978569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 106078569bb4SMatthew G. Knepley 106178569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 106278569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 106378569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 106478569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 106578569bb4SMatthew G. Knepley if (bs == 1) { 10669566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vNatural, &varray)); 1067a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]); 10689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vNatural, &varray)); 106978569bb4SMatthew G. Knepley } else { 107078569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 10719566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs)); 10729566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 10739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(vComp, &varray)); 1074a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]); 10759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(vComp, &varray)); 10769566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 107778569bb4SMatthew G. Knepley } 107878569bb4SMatthew G. Knepley } 107978569bb4SMatthew G. Knepley csxs += csSize[set]; 108078569bb4SMatthew G. Knepley } 10819566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 10829566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 10839566063dSJacob Faibussowitsch if (useNatural) PetscCall(DMRestoreGlobalVector(dm,&vNatural)); 108478569bb4SMatthew G. Knepley PetscFunctionReturn(0); 108578569bb4SMatthew G. Knepley } 108678569bb4SMatthew G. Knepley 108778569bb4SMatthew G. Knepley /* 108878569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 108978569bb4SMatthew G. Knepley 109078569bb4SMatthew G. Knepley Collective on v 109178569bb4SMatthew G. Knepley 109278569bb4SMatthew G. Knepley Input Parameters: 109378569bb4SMatthew G. Knepley + v - The vector to be written 109478569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 10956823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1) 10966823f3c5SBlaise Bourdin - offset - the location of the variable in the file 109778569bb4SMatthew G. Knepley 109878569bb4SMatthew G. Knepley Notes: 109978569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 110078569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 110178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 110278569bb4SMatthew G. Knepley amongst all zonal variables. 110378569bb4SMatthew G. Knepley 110478569bb4SMatthew G. Knepley Level: beginner 110578569bb4SMatthew G. Knepley 1106db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()` 110778569bb4SMatthew G. Knepley */ 11086823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) 110978569bb4SMatthew G. Knepley { 111078569bb4SMatthew G. Knepley MPI_Comm comm; 111178569bb4SMatthew G. Knepley PetscMPIInt size; 111278569bb4SMatthew G. Knepley DM dm; 111378569bb4SMatthew G. Knepley Vec vNatural, vComp; 111422a7544eSVaclav Hapla PetscScalar *varray; 111578569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 111678569bb4SMatthew G. Knepley PetscBool useNatural; 111778569bb4SMatthew G. Knepley IS compIS; 111878569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 111978569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 112078569bb4SMatthew G. Knepley 112178569bb4SMatthew G. Knepley PetscFunctionBegin; 11229566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)v,&comm)); 11239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 11249566063dSJacob Faibussowitsch PetscCall(VecGetDM(v, &dm)); 11259566063dSJacob Faibussowitsch PetscCall(DMGetUseNatural(dm, &useNatural)); 112678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 11279566063dSJacob Faibussowitsch if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural)); 112878569bb4SMatthew G. Knepley else {vNatural = v;} 11296823f3c5SBlaise Bourdin 113078569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 113178569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 113278569bb4SMatthew G. Knepley We assume that they are stored sequentially 113378569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 1134a5b23f4aSJose E. Roman but once the vector has been reordered to natural size, we cannot use the label information 113578569bb4SMatthew G. Knepley to figure out what to save where. */ 113678569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 11379566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize)); 1138a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID); 113978569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 114078569bb4SMatthew G. Knepley ex_block block; 114178569bb4SMatthew G. Knepley 114278569bb4SMatthew G. Knepley block.id = csID[set]; 114378569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 1144a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_block_param,exoid, &block); 114578569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 114678569bb4SMatthew G. Knepley } 11479566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe)); 11489566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(vNatural, &bs)); 11499566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS)); 115078569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 115178569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 115278569bb4SMatthew G. Knepley 115378569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 115478569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 115578569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 115678569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 115778569bb4SMatthew G. Knepley if (bs == 1) { 11589566063dSJacob Faibussowitsch PetscCall(VecGetArray(vNatural, &varray)); 1159a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]); 11609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vNatural, &varray)); 116178569bb4SMatthew G. Knepley } else { 116278569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 11639566063dSJacob Faibussowitsch PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs)); 11649566063dSJacob Faibussowitsch PetscCall(VecGetSubVector(vNatural, compIS, &vComp)); 11659566063dSJacob Faibussowitsch PetscCall(VecGetArray(vComp, &varray)); 1166a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]); 11679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vComp, &varray)); 11689566063dSJacob Faibussowitsch PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp)); 116978569bb4SMatthew G. Knepley } 117078569bb4SMatthew G. Knepley } 117178569bb4SMatthew G. Knepley csxs += csSize[set]; 117278569bb4SMatthew G. Knepley } 11739566063dSJacob Faibussowitsch PetscCall(PetscFree2(csID, csSize)); 11749566063dSJacob Faibussowitsch if (bs > 1) PetscCall(ISDestroy(&compIS)); 117578569bb4SMatthew G. Knepley if (useNatural) { 11769566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v)); 11779566063dSJacob Faibussowitsch PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v)); 11789566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &vNatural)); 117978569bb4SMatthew G. Knepley } 118078569bb4SMatthew G. Knepley PetscFunctionReturn(0); 118178569bb4SMatthew G. Knepley } 1182b53c8456SSatish Balay #endif 118378569bb4SMatthew G. Knepley 11841e50132fSMatthew G. Knepley /*@ 11851e50132fSMatthew G. Knepley PetscViewerExodusIIGetId - Get the file id of the ExodusII file 11861e50132fSMatthew G. Knepley 11871e50132fSMatthew G. Knepley Logically Collective on PetscViewer 11881e50132fSMatthew G. Knepley 11891e50132fSMatthew G. Knepley Input Parameter: 11901e50132fSMatthew G. Knepley . viewer - the PetscViewer 11911e50132fSMatthew G. Knepley 11921e50132fSMatthew G. Knepley Output Parameter: 1193d8d19677SJose E. Roman . exoid - The ExodusII file id 11941e50132fSMatthew G. Knepley 11951e50132fSMatthew G. Knepley Level: intermediate 11961e50132fSMatthew G. Knepley 1197db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 11981e50132fSMatthew G. Knepley @*/ 11991e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) 12001e50132fSMatthew G. Knepley { 12011e50132fSMatthew G. Knepley PetscFunctionBegin; 12021e50132fSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1203cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetId_C",(PetscViewer,int*),(viewer,exoid)); 12046823f3c5SBlaise Bourdin PetscFunctionReturn(0); 12056823f3c5SBlaise Bourdin } 12066823f3c5SBlaise Bourdin 12076823f3c5SBlaise Bourdin /*@ 12086823f3c5SBlaise Bourdin PetscViewerExodusIISetOrder - Set the elements order in the exodusII file. 12096823f3c5SBlaise Bourdin 12106823f3c5SBlaise Bourdin Collective 12116823f3c5SBlaise Bourdin 12126823f3c5SBlaise Bourdin Input Parameters: 12136823f3c5SBlaise Bourdin + viewer - the viewer 121498a6a78aSPierre Jolivet - order - elements order 12156823f3c5SBlaise Bourdin 12166823f3c5SBlaise Bourdin Output Parameter: 12176823f3c5SBlaise Bourdin 12186823f3c5SBlaise Bourdin Level: beginner 12196823f3c5SBlaise Bourdin 12206823f3c5SBlaise Bourdin Note: 12216823f3c5SBlaise Bourdin 1222db781477SPatrick Sanan .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 12236823f3c5SBlaise Bourdin @*/ 12246823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order) 12256823f3c5SBlaise Bourdin { 12266823f3c5SBlaise Bourdin PetscFunctionBegin; 12276823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1228cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerSetOrder_C",(PetscViewer,PetscInt),(viewer,order)); 12296823f3c5SBlaise Bourdin PetscFunctionReturn(0); 12306823f3c5SBlaise Bourdin } 12316823f3c5SBlaise Bourdin 12326823f3c5SBlaise Bourdin /*@ 12336823f3c5SBlaise Bourdin PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file. 12346823f3c5SBlaise Bourdin 12356823f3c5SBlaise Bourdin Collective 12366823f3c5SBlaise Bourdin 12376823f3c5SBlaise Bourdin Input Parameters: 12386823f3c5SBlaise Bourdin + viewer - the viewer 123998a6a78aSPierre Jolivet - order - elements order 12406823f3c5SBlaise Bourdin 12416823f3c5SBlaise Bourdin Output Parameter: 12426823f3c5SBlaise Bourdin 12436823f3c5SBlaise Bourdin Level: beginner 12446823f3c5SBlaise Bourdin 12456823f3c5SBlaise Bourdin Note: 12466823f3c5SBlaise Bourdin 1247db781477SPatrick Sanan .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()` 12486823f3c5SBlaise Bourdin @*/ 12496823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order) 12506823f3c5SBlaise Bourdin { 12516823f3c5SBlaise Bourdin PetscFunctionBegin; 12526823f3c5SBlaise Bourdin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1253cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerGetOrder_C",(PetscViewer,PetscInt*),(viewer,order)); 12546823f3c5SBlaise Bourdin PetscFunctionReturn(0); 12556823f3c5SBlaise Bourdin } 12566823f3c5SBlaise Bourdin 12576823f3c5SBlaise Bourdin /*@C 12586823f3c5SBlaise Bourdin PetscViewerExodusIIOpen - Opens a file for ExodusII input/output. 12596823f3c5SBlaise Bourdin 12606823f3c5SBlaise Bourdin Collective 12616823f3c5SBlaise Bourdin 12626823f3c5SBlaise Bourdin Input Parameters: 12636823f3c5SBlaise Bourdin + comm - MPI communicator 12646823f3c5SBlaise Bourdin . name - name of file 12656823f3c5SBlaise Bourdin - type - type of file 12666823f3c5SBlaise Bourdin $ FILE_MODE_WRITE - create new file for binary output 12676823f3c5SBlaise Bourdin $ FILE_MODE_READ - open existing file for binary input 12686823f3c5SBlaise Bourdin $ FILE_MODE_APPEND - open existing file for binary output 12696823f3c5SBlaise Bourdin 12706823f3c5SBlaise Bourdin Output Parameter: 12716823f3c5SBlaise Bourdin . exo - PetscViewer for Exodus II input/output to use with the specified file 12726823f3c5SBlaise Bourdin 12736823f3c5SBlaise Bourdin Level: beginner 12746823f3c5SBlaise Bourdin 12756823f3c5SBlaise Bourdin Note: 12766823f3c5SBlaise Bourdin This PetscViewer should be destroyed with PetscViewerDestroy(). 12776823f3c5SBlaise Bourdin 1278db781477SPatrick Sanan .seealso: `PetscViewerPushFormat()`, `PetscViewerDestroy()`, 1279db781477SPatrick Sanan `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 12806823f3c5SBlaise Bourdin @*/ 12816823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) 12826823f3c5SBlaise Bourdin { 12836823f3c5SBlaise Bourdin PetscFunctionBegin; 12849566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, exo)); 12859566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII)); 12869566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*exo, type)); 12879566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*exo, name)); 12889566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*exo)); 12891e50132fSMatthew G. Knepley PetscFunctionReturn(0); 12901e50132fSMatthew G. Knepley } 12911e50132fSMatthew G. Knepley 129233751fbdSMichael Lange /*@C 1293eaf898f9SPatrick Sanan DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 129433751fbdSMichael Lange 1295d083f849SBarry Smith Collective 129633751fbdSMichael Lange 129733751fbdSMichael Lange Input Parameters: 129833751fbdSMichael Lange + comm - The MPI communicator 129933751fbdSMichael Lange . filename - The name of the ExodusII file 130033751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 130133751fbdSMichael Lange 130233751fbdSMichael Lange Output Parameter: 130333751fbdSMichael Lange . dm - The DM object representing the mesh 130433751fbdSMichael Lange 130533751fbdSMichael Lange Level: beginner 130633751fbdSMichael Lange 1307db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()` 130833751fbdSMichael Lange @*/ 130933751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 131033751fbdSMichael Lange { 131133751fbdSMichael Lange PetscMPIInt rank; 131233751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1313e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 131433751fbdSMichael Lange float version; 131533751fbdSMichael Lange #endif 131633751fbdSMichael Lange 131733751fbdSMichael Lange PetscFunctionBegin; 131833751fbdSMichael Lange PetscValidCharPointer(filename, 2); 13199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 132033751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1321dd400576SPatrick Sanan if (rank == 0) { 132233751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 132308401ef6SPierre Jolivet PetscCheck(exoid > 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 132433751fbdSMichael Lange } 13259566063dSJacob Faibussowitsch PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm)); 1326a74df02fSJacob Faibussowitsch if (rank == 0) {PetscStackCallStandard(ex_close,exoid);} 1327b458e8f1SJose E. Roman PetscFunctionReturn(0); 132833751fbdSMichael Lange #else 132933751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 133033751fbdSMichael Lange #endif 133133751fbdSMichael Lange } 133233751fbdSMichael Lange 13338f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII) 13346823f3c5SBlaise Bourdin static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct) 13358f861fbcSMatthew G. Knepley { 13368f861fbcSMatthew G. Knepley PetscBool flg; 13378f861fbcSMatthew G. Knepley 13388f861fbcSMatthew G. Knepley PetscFunctionBegin; 13398f861fbcSMatthew G. Knepley *ct = DM_POLYTOPE_UNKNOWN; 13409566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI", &flg)); 13418f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;} 13429566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TRI3", &flg)); 13438f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;} 13449566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD", &flg)); 13458f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;} 13469566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg)); 13478f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;} 13489566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg)); 13498f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;} 13509566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TETRA", &flg)); 13518f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;} 13529566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "TET4", &flg)); 13538f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;} 13549566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg)); 13558f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_TRI_PRISM; goto done;} 13569566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX", &flg)); 13578f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;} 13589566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEX8", &flg)); 13598f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;} 13609566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg)); 13618f861fbcSMatthew G. Knepley if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;} 136298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type); 13638f861fbcSMatthew G. Knepley done: 13648f861fbcSMatthew G. Knepley PetscFunctionReturn(0); 13658f861fbcSMatthew G. Knepley } 13668f861fbcSMatthew G. Knepley #endif 13678f861fbcSMatthew G. Knepley 1368552f7358SJed Brown /*@ 136933751fbdSMichael Lange DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 1370552f7358SJed Brown 1371d083f849SBarry Smith Collective 1372552f7358SJed Brown 1373552f7358SJed Brown Input Parameters: 1374552f7358SJed Brown + comm - The MPI communicator 1375552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 1376552f7358SJed Brown - interpolate - Create faces and edges in the mesh 1377552f7358SJed Brown 1378552f7358SJed Brown Output Parameter: 1379552f7358SJed Brown . dm - The DM object representing the mesh 1380552f7358SJed Brown 1381552f7358SJed Brown Level: beginner 1382552f7358SJed Brown 1383db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()` 1384552f7358SJed Brown @*/ 1385552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 1386552f7358SJed Brown { 1387552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 1388552f7358SJed Brown PetscMPIInt num_proc, rank; 1389ae9eebabSLisandro Dalcin DMLabel cellSets = NULL, faceSets = NULL, vertSets = NULL; 1390552f7358SJed Brown PetscSection coordSection; 1391552f7358SJed Brown Vec coordinates; 1392552f7358SJed Brown PetscScalar *coords; 1393552f7358SJed Brown PetscInt coordSize, v; 1394552f7358SJed Brown /* Read from ex_get_init() */ 1395552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 13965f80ce2aSJacob Faibussowitsch int dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0; 1397552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 1398552f7358SJed Brown #endif 1399552f7358SJed Brown 1400552f7358SJed Brown PetscFunctionBegin; 1401552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 14029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 14039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &num_proc)); 14049566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 14059566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMPLEX)); 1406a5b23f4aSJose E. Roman /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */ 1407dd400576SPatrick Sanan if (rank == 0) { 14089566063dSJacob Faibussowitsch PetscCall(PetscMemzero(title,PETSC_MAX_PATH_LEN+1)); 1409a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_init,exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs); 141028b400f6SJacob Faibussowitsch PetscCheck(num_cs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set"); 1411552f7358SJed Brown } 14129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm)); 14139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm)); 14149566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) *dm, title)); 14159566063dSJacob Faibussowitsch PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices)); 1416412e9a14SMatthew G. Knepley /* We do not want this label automatically computed, instead we compute it here */ 14179566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(*dm, "celltype")); 1418552f7358SJed Brown 1419552f7358SJed Brown /* Read cell sets information */ 1420dd400576SPatrick Sanan if (rank == 0) { 1421552f7358SJed Brown PetscInt *cone; 1422e8f6893fSMatthew G. Knepley int c, cs, ncs, c_loc, v, v_loc; 1423552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 1424e8f6893fSMatthew G. Knepley int *cs_id, *cs_order; 1425552f7358SJed Brown /* Read from ex_get_elem_block() */ 1426552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 1427e8f6893fSMatthew G. Knepley int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr; 1428552f7358SJed Brown /* Read from ex_get_elem_conn() */ 1429552f7358SJed Brown int *cs_connect; 1430552f7358SJed Brown 1431552f7358SJed Brown /* Get cell sets IDs */ 14329566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order)); 1433a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, cs_id); 1434552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 1435552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 1436e8f6893fSMatthew G. Knepley /* Check for a hybrid mesh */ 1437e8f6893fSMatthew G. Knepley for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) { 14388f861fbcSMatthew G. Knepley DMPolytopeType ct; 14398f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 14408f861fbcSMatthew G. Knepley 14419566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1442a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type); 14439566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 14448f861fbcSMatthew G. Knepley dim = PetscMax(dim, DMPolytopeTypeGetDim(ct)); 1445a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr); 14468f861fbcSMatthew G. Knepley switch (ct) { 14478f861fbcSMatthew G. Knepley case DM_POLYTOPE_TRI_PRISM: 1448e8f6893fSMatthew G. Knepley cs_order[cs] = cs; 1449e8f6893fSMatthew G. Knepley ++num_hybrid; 1450e8f6893fSMatthew G. Knepley break; 1451e8f6893fSMatthew G. Knepley default: 1452e8f6893fSMatthew G. Knepley for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1]; 1453e8f6893fSMatthew G. Knepley cs_order[cs-num_hybrid] = cs; 1454e8f6893fSMatthew G. Knepley } 1455e8f6893fSMatthew G. Knepley } 1456552f7358SJed Brown /* First set sizes */ 1457e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 14588f861fbcSMatthew G. Knepley DMPolytopeType ct; 14598f861fbcSMatthew G. Knepley char elem_type[PETSC_MAX_PATH_LEN]; 1460e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 14618f861fbcSMatthew G. Knepley 14629566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(elem_type, sizeof(elem_type))); 1463a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type); 14649566063dSJacob Faibussowitsch PetscCall(ExodusGetCellType_Internal(elem_type, &ct)); 1465a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr); 1466552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 14679566063dSJacob Faibussowitsch PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell)); 14689566063dSJacob Faibussowitsch PetscCall(DMPlexSetCellType(*dm, c, ct)); 1469552f7358SJed Brown } 1470552f7358SJed Brown } 14719566063dSJacob Faibussowitsch for (v = numCells; v < numCells+numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT)); 14729566063dSJacob Faibussowitsch PetscCall(DMSetUp(*dm)); 1473e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1474e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 1475a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr); 14769566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone)); 1477a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_conn,exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL); 1478eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 1479552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1480636e64ffSStefano Zampini DMPolytopeType ct; 1481636e64ffSStefano Zampini 1482552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 1483552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 1484552f7358SJed Brown } 14859566063dSJacob Faibussowitsch PetscCall(DMPlexGetCellType(*dm, c, &ct)); 14869566063dSJacob Faibussowitsch PetscCall(DMPlexInvertCell(ct, cone)); 14879566063dSJacob Faibussowitsch PetscCall(DMPlexSetCone(*dm, c, cone)); 14889566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs])); 1489552f7358SJed Brown } 14909566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_connect,cone)); 1491552f7358SJed Brown } 14929566063dSJacob Faibussowitsch PetscCall(PetscFree2(cs_id, cs_order)); 1493552f7358SJed Brown } 14948f861fbcSMatthew G. Knepley { 14958f861fbcSMatthew G. Knepley PetscInt ints[] = {dim, dimEmbed}; 14968f861fbcSMatthew G. Knepley 14979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm)); 14989566063dSJacob Faibussowitsch PetscCall(DMSetDimension(*dm, ints[0])); 14999566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(*dm, ints[1])); 15008f861fbcSMatthew G. Knepley dim = ints[0]; 15018f861fbcSMatthew G. Knepley dimEmbed = ints[1]; 15028f861fbcSMatthew G. Knepley } 15039566063dSJacob Faibussowitsch PetscCall(DMPlexSymmetrize(*dm)); 15049566063dSJacob Faibussowitsch PetscCall(DMPlexStratify(*dm)); 1505552f7358SJed Brown if (interpolate) { 15065fd9971aSMatthew G. Knepley DM idm; 1507552f7358SJed Brown 15089566063dSJacob Faibussowitsch PetscCall(DMPlexInterpolate(*dm, &idm)); 15099566063dSJacob Faibussowitsch PetscCall(DMDestroy(dm)); 1510552f7358SJed Brown *dm = idm; 1511552f7358SJed Brown } 1512552f7358SJed Brown 1513552f7358SJed Brown /* Create vertex set label */ 1514dd400576SPatrick Sanan if (rank == 0 && (num_vs > 0)) { 1515552f7358SJed Brown int vs, v; 1516552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 1517552f7358SJed Brown int *vs_id; 1518552f7358SJed Brown /* Read from ex_get_node_set_param() */ 1519f958083aSBlaise Bourdin int num_vertex_in_set; 1520552f7358SJed Brown /* Read from ex_get_node_set() */ 1521552f7358SJed Brown int *vs_vertex_list; 1522552f7358SJed Brown 1523552f7358SJed Brown /* Get vertex set ids */ 15249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vs, &vs_id)); 1525a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_ids,exoid, EX_NODE_SET, vs_id); 1526552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 1527a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_set_param,exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL); 15289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list)); 1529a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_set,exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL); 1530552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 15319566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs])); 1532552f7358SJed Brown } 15339566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_vertex_list)); 1534552f7358SJed Brown } 15359566063dSJacob Faibussowitsch PetscCall(PetscFree(vs_id)); 1536552f7358SJed Brown } 1537552f7358SJed Brown /* Read coordinates */ 15389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateSection(*dm, &coordSection)); 15399566063dSJacob Faibussowitsch PetscCall(PetscSectionSetNumFields(coordSection, 1)); 15409566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed)); 15419566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices)); 1542552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 15439566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed)); 15449566063dSJacob Faibussowitsch PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed)); 1545552f7358SJed Brown } 15469566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(coordSection)); 15479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize)); 15489566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates)); 15499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates")); 15509566063dSJacob Faibussowitsch PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE)); 15519566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(coordinates, dimEmbed)); 15529566063dSJacob Faibussowitsch PetscCall(VecSetType(coordinates,VECSTANDARD)); 15539566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 1554dd400576SPatrick Sanan if (rank == 0) { 1555e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 1556552f7358SJed Brown 15579566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z)); 1558a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_coord,exoid, x, y, z); 15598f861fbcSMatthew G. Knepley if (dimEmbed > 0) { 15608f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+0] = x[v]; 15610d644c17SKarl Rupp } 15628f861fbcSMatthew G. Knepley if (dimEmbed > 1) { 15638f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+1] = y[v]; 15640d644c17SKarl Rupp } 15658f861fbcSMatthew G. Knepley if (dimEmbed > 2) { 15668f861fbcSMatthew G. Knepley for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+2] = z[v]; 15670d644c17SKarl Rupp } 15689566063dSJacob Faibussowitsch PetscCall(PetscFree3(x,y,z)); 1569552f7358SJed Brown } 15709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 15719566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(*dm, coordinates)); 15729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&coordinates)); 1573552f7358SJed Brown 1574552f7358SJed Brown /* Create side set label */ 1575dd400576SPatrick Sanan if (rank == 0 && interpolate && (num_fs > 0)) { 1576552f7358SJed Brown int fs, f, voff; 1577552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 1578552f7358SJed Brown int *fs_id; 1579552f7358SJed Brown /* Read from ex_get_side_set_param() */ 1580f958083aSBlaise Bourdin int num_side_in_set; 1581552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 1582552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 1583ef073283Smichael_afanasiev /* Read side set labels */ 1584ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH+1]; 1585552f7358SJed Brown 1586552f7358SJed Brown /* Get side set ids */ 15879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(num_fs, &fs_id)); 1588a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_ids,exoid, EX_SIDE_SET, fs_id); 1589552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 1590a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_set_param,exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL); 15919566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list)); 1592a74df02fSJacob Faibussowitsch PetscStackCallStandard(ex_get_side_set_node_list,exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list); 1593ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 1594ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1595552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 15960298fd71SBarry Smith const PetscInt *faces = NULL; 1597552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1598552f7358SJed Brown PetscInt faceVertices[4], v; 1599552f7358SJed Brown 160008401ef6SPierre Jolivet PetscCheck(faceSize <= 4,comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 1601552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 1602552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 1603552f7358SJed Brown } 16049566063dSJacob Faibussowitsch PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 160508401ef6SPierre Jolivet PetscCheck(numFaces == 1,comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 16069566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs])); 1607ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 1608ef073283Smichael_afanasiev if (!fs_name_err) { 16099566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs])); 1610ef073283Smichael_afanasiev } 16119566063dSJacob Faibussowitsch PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces)); 1612552f7358SJed Brown } 16139566063dSJacob Faibussowitsch PetscCall(PetscFree2(fs_vertex_count_list,fs_vertex_list)); 1614552f7358SJed Brown } 16159566063dSJacob Faibussowitsch PetscCall(PetscFree(fs_id)); 1616552f7358SJed Brown } 1617ae9eebabSLisandro Dalcin 1618ae9eebabSLisandro Dalcin { /* Create Cell/Face/Vertex Sets labels at all processes */ 1619ae9eebabSLisandro Dalcin enum {n = 3}; 1620ae9eebabSLisandro Dalcin PetscBool flag[n]; 1621ae9eebabSLisandro Dalcin 1622ae9eebabSLisandro Dalcin flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE; 1623ae9eebabSLisandro Dalcin flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE; 1624ae9eebabSLisandro Dalcin flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE; 16259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm)); 16269566063dSJacob Faibussowitsch if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets")); 16279566063dSJacob Faibussowitsch if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets")); 16289566063dSJacob Faibussowitsch if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets")); 1629ae9eebabSLisandro Dalcin } 1630b458e8f1SJose E. Roman PetscFunctionReturn(0); 1631552f7358SJed Brown #else 1632552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1633552f7358SJed Brown #endif 1634552f7358SJed Brown } 1635