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 3000373969SVaclav Hapla .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 { 471e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 481e50132fSMatthew G. Knepley PetscErrorCode ierr; 491e50132fSMatthew G. Knepley 501e50132fSMatthew G. Knepley PetscFunctionBegin; 511e50132fSMatthew G. Knepley if (exo->filename) {ierr = PetscViewerASCIIPrintf(viewer, "Filename: %s\n", exo->filename);CHKERRQ(ierr);} 521e50132fSMatthew G. Knepley PetscFunctionReturn(0); 531e50132fSMatthew G. Knepley } 541e50132fSMatthew G. Knepley 551e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v) 561e50132fSMatthew G. Knepley { 571e50132fSMatthew G. Knepley PetscErrorCode ierr; 581e50132fSMatthew G. Knepley 591e50132fSMatthew G. Knepley PetscFunctionBegin; 601e50132fSMatthew G. Knepley ierr = PetscOptionsHead(PetscOptionsObject, "ExodusII PetscViewer Options");CHKERRQ(ierr); 611e50132fSMatthew G. Knepley ierr = PetscOptionsTail();CHKERRQ(ierr); 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 PetscErrorCode ierr; 751e50132fSMatthew G. Knepley 761e50132fSMatthew G. Knepley PetscFunctionBegin; 771e50132fSMatthew G. Knepley if (exo->exoid >= 0) {ierr = ex_close(exo->exoid);CHKERRQ(ierr);} 781e50132fSMatthew G. Knepley ierr = PetscFree(exo);CHKERRQ(ierr); 791e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 801e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 811e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 821e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerExodusIIGetId",NULL);CHKERRQ(ierr); 831e50132fSMatthew G. Knepley PetscFunctionReturn(0); 841e50132fSMatthew G. Knepley } 851e50132fSMatthew G. Knepley 861e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[]) 871e50132fSMatthew G. Knepley { 881e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 891e50132fSMatthew G. Knepley PetscMPIInt rank; 901e50132fSMatthew G. Knepley int CPU_word_size, IO_word_size, EXO_mode; 911e50132fSMatthew G. Knepley PetscErrorCode ierr; 921e50132fSMatthew G. Knepley 931e50132fSMatthew G. Knepley 941e50132fSMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank);CHKERRQ(ierr); 951e50132fSMatthew G. Knepley CPU_word_size = sizeof(PetscReal); 961e50132fSMatthew G. Knepley IO_word_size = sizeof(PetscReal); 971e50132fSMatthew G. Knepley 981e50132fSMatthew G. Knepley PetscFunctionBegin; 991e50132fSMatthew G. Knepley if (exo->exoid >= 0) {ex_close(exo->exoid); exo->exoid = -1;} 1001e50132fSMatthew G. Knepley if (exo->filename) {ierr = PetscFree(exo->filename);CHKERRQ(ierr);} 1011e50132fSMatthew G. Knepley ierr = PetscStrallocpy(name, &exo->filename);CHKERRQ(ierr); 1021e50132fSMatthew G. Knepley /* Create or open the file collectively */ 1031e50132fSMatthew G. Knepley switch (exo->btype) { 1041e50132fSMatthew G. Knepley case FILE_MODE_READ: 1051e50132fSMatthew G. Knepley EXO_mode = EX_CLOBBER; 1061e50132fSMatthew G. Knepley break; 1071e50132fSMatthew G. Knepley case FILE_MODE_APPEND: 1081e50132fSMatthew G. Knepley EXO_mode = EX_CLOBBER; 1091e50132fSMatthew G. Knepley break; 1101e50132fSMatthew G. Knepley case FILE_MODE_WRITE: 1111e50132fSMatthew G. Knepley EXO_mode = EX_CLOBBER; 1121e50132fSMatthew G. Knepley break; 1131e50132fSMatthew G. Knepley default: 1141e50132fSMatthew G. Knepley SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 1151e50132fSMatthew G. Knepley } 1161e50132fSMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 1171e50132fSMatthew G. Knepley EXO_mode += EX_ALL_INT64_API; 1181e50132fSMatthew G. Knepley #endif 1191e50132fSMatthew G. Knepley exo->exoid = ex_create(name, EXO_mode, &CPU_word_size, &IO_word_size); 1201e50132fSMatthew G. Knepley if (exo->exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", name); 1211e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1221e50132fSMatthew G. Knepley } 1231e50132fSMatthew G. Knepley 1241e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name) 1251e50132fSMatthew G. Knepley { 1261e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1271e50132fSMatthew G. Knepley 1281e50132fSMatthew G. Knepley PetscFunctionBegin; 1291e50132fSMatthew G. Knepley *name = exo->filename; 1301e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1311e50132fSMatthew G. Knepley } 1321e50132fSMatthew G. Knepley 1331e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type) 1341e50132fSMatthew G. Knepley { 1351e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1361e50132fSMatthew G. Knepley 1371e50132fSMatthew G. Knepley PetscFunctionBegin; 1381e50132fSMatthew G. Knepley exo->btype = type; 1391e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1401e50132fSMatthew G. Knepley } 1411e50132fSMatthew G. Knepley 1421e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type) 1431e50132fSMatthew G. Knepley { 1441e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1451e50132fSMatthew G. Knepley 1461e50132fSMatthew G. Knepley PetscFunctionBegin; 1471e50132fSMatthew G. Knepley *type = exo->btype; 1481e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1491e50132fSMatthew G. Knepley } 1501e50132fSMatthew G. Knepley 1511e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid) 1521e50132fSMatthew G. Knepley { 1531e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data; 1541e50132fSMatthew G. Knepley 1551e50132fSMatthew G. Knepley PetscFunctionBegin; 1561e50132fSMatthew G. Knepley *exoid = exo->exoid; 1571e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1581e50132fSMatthew G. Knepley } 1591e50132fSMatthew G. Knepley 1601e50132fSMatthew G. Knepley /*@C 1611e50132fSMatthew G. Knepley PetscViewerExodusIIOpen - Opens a file for ExodusII input/output. 1621e50132fSMatthew G. Knepley 1631e50132fSMatthew G. Knepley Collective 1641e50132fSMatthew G. Knepley 1651e50132fSMatthew G. Knepley Input Parameters: 1661e50132fSMatthew G. Knepley + comm - MPI communicator 1671e50132fSMatthew G. Knepley . name - name of file 1681e50132fSMatthew G. Knepley - type - type of file 1691e50132fSMatthew G. Knepley $ FILE_MODE_WRITE - create new file for binary output 1701e50132fSMatthew G. Knepley $ FILE_MODE_READ - open existing file for binary input 1711e50132fSMatthew G. Knepley $ FILE_MODE_APPEND - open existing file for binary output 1721e50132fSMatthew G. Knepley 1731e50132fSMatthew G. Knepley Output Parameter: 17400373969SVaclav Hapla . exo - PetscViewer for Exodus II input/output to use with the specified file 1751e50132fSMatthew G. Knepley 1761e50132fSMatthew G. Knepley Level: beginner 1771e50132fSMatthew G. Knepley 1781e50132fSMatthew G. Knepley Note: 1791e50132fSMatthew G. Knepley This PetscViewer should be destroyed with PetscViewerDestroy(). 1801e50132fSMatthew G. Knepley 1811e50132fSMatthew G. Knepley 18200373969SVaclav Hapla .seealso: PetscViewerPushFormat(), PetscViewerDestroy(), 18300373969SVaclav Hapla DMLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 1841e50132fSMatthew G. Knepley @*/ 1851e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) 1861e50132fSMatthew G. Knepley { 1871e50132fSMatthew G. Knepley PetscErrorCode ierr; 1881e50132fSMatthew G. Knepley 1891e50132fSMatthew G. Knepley PetscFunctionBegin; 1901e50132fSMatthew G. Knepley ierr = PetscViewerCreate(comm, exo);CHKERRQ(ierr); 1911e50132fSMatthew G. Knepley ierr = PetscViewerSetType(*exo, PETSCVIEWEREXODUSII);CHKERRQ(ierr); 1921e50132fSMatthew G. Knepley ierr = PetscViewerFileSetMode(*exo, type);CHKERRQ(ierr); 1931e50132fSMatthew G. Knepley ierr = PetscViewerFileSetName(*exo, name);CHKERRQ(ierr); 1941e50132fSMatthew G. Knepley ierr = PetscViewerSetFromOptions(*exo);CHKERRQ(ierr); 1951e50132fSMatthew G. Knepley PetscFunctionReturn(0); 1961e50132fSMatthew G. Knepley } 1971e50132fSMatthew G. Knepley 1981e50132fSMatthew G. Knepley /*MC 19900373969SVaclav Hapla PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file 2001e50132fSMatthew G. Knepley 2011e50132fSMatthew G. Knepley 20200373969SVaclav Hapla .seealso: PetscViewerExodusIIOpen(), PetscViewerCreate(), PETSCVIEWERBINARY, PETSCVIEWERHDF5, DMView(), 2031e50132fSMatthew G. Knepley PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 2041e50132fSMatthew G. Knepley 2051e50132fSMatthew G. Knepley Level: beginner 2061e50132fSMatthew G. Knepley M*/ 2071e50132fSMatthew G. Knepley 2081e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v) 2091e50132fSMatthew G. Knepley { 2101e50132fSMatthew G. Knepley PetscViewer_ExodusII *exo; 2111e50132fSMatthew G. Knepley PetscErrorCode ierr; 2121e50132fSMatthew G. Knepley 2131e50132fSMatthew G. Knepley PetscFunctionBegin; 2141e50132fSMatthew G. Knepley ierr = PetscNewLog(v,&exo);CHKERRQ(ierr); 2151e50132fSMatthew G. Knepley 2161e50132fSMatthew G. Knepley v->data = (void*) exo; 2171e50132fSMatthew G. Knepley v->ops->destroy = PetscViewerDestroy_ExodusII; 2181e50132fSMatthew G. Knepley v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII; 2191e50132fSMatthew G. Knepley v->ops->setup = PetscViewerSetUp_ExodusII; 2201e50132fSMatthew G. Knepley v->ops->view = PetscViewerView_ExodusII; 2211e50132fSMatthew G. Knepley v->ops->flush = 0; 2221e50132fSMatthew G. Knepley exo->btype = (PetscFileMode) -1; 2231e50132fSMatthew G. Knepley exo->filename = 0; 2241e50132fSMatthew G. Knepley exo->exoid = -1; 2251e50132fSMatthew G. Knepley 2261e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII);CHKERRQ(ierr); 2271e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII);CHKERRQ(ierr); 2281e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII);CHKERRQ(ierr); 2291e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII);CHKERRQ(ierr); 2301e50132fSMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerExodusIIGetId_C",PetscViewerExodusIIGetId_ExodusII);CHKERRQ(ierr); 2311e50132fSMatthew G. Knepley PetscFunctionReturn(0); 2321e50132fSMatthew G. Knepley } 2331e50132fSMatthew G. Knepley 2341e50132fSMatthew G. Knepley /* 23578569bb4SMatthew G. Knepley EXOGetVarIndex - Locate a result in an exodus file based on its name 23678569bb4SMatthew G. Knepley 23778569bb4SMatthew G. Knepley Not collective 23878569bb4SMatthew G. Knepley 23978569bb4SMatthew G. Knepley Input Parameters: 24078569bb4SMatthew G. Knepley + exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 24178569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK 24278569bb4SMatthew G. Knepley - name - the name of the result 24378569bb4SMatthew G. Knepley 24478569bb4SMatthew G. Knepley Output Parameters: 24578569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result 24678569bb4SMatthew G. Knepley 24778569bb4SMatthew G. Knepley Notes: 24878569bb4SMatthew G. Knepley The exodus variable index is obtained by comparing name and the 24978569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance if name is "V" 25078569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 25178569bb4SMatthew G. Knepley amongst all variables of type obj_type. 25278569bb4SMatthew G. Knepley 25378569bb4SMatthew G. Knepley Level: beginner 25478569bb4SMatthew G. Knepley 255cd921712SStefano Zampini .seealso: EXOGetVarIndex(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO() 25678569bb4SMatthew G. Knepley */ 25778569bb4SMatthew G. Knepley static PetscErrorCode EXOGetVarIndex_Private(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) 25878569bb4SMatthew G. Knepley { 2596de834b4SMatthew G. Knepley int num_vars, i, j; 26078569bb4SMatthew G. Knepley char ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1]; 26178569bb4SMatthew G. Knepley const int num_suffix = 5; 26278569bb4SMatthew G. Knepley char *suffix[5]; 26378569bb4SMatthew G. Knepley PetscBool flg; 26478569bb4SMatthew G. Knepley PetscErrorCode ierr; 26578569bb4SMatthew G. Knepley 26678569bb4SMatthew G. Knepley PetscFunctionBegin; 26778569bb4SMatthew G. Knepley suffix[0] = (char *) ""; 26878569bb4SMatthew G. Knepley suffix[1] = (char *) "_X"; 26978569bb4SMatthew G. Knepley suffix[2] = (char *) "_XX"; 27078569bb4SMatthew G. Knepley suffix[3] = (char *) "_1"; 27178569bb4SMatthew G. Knepley suffix[4] = (char *) "_11"; 27278569bb4SMatthew G. Knepley 27378569bb4SMatthew G. Knepley *varIndex = 0; 2746de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_variable_param,(exoid, obj_type, &num_vars)); 27578569bb4SMatthew G. Knepley for (i = 0; i < num_vars; ++i) { 2766de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_variable_name,(exoid, obj_type, i+1, var_name)); 27778569bb4SMatthew G. Knepley for (j = 0; j < num_suffix; ++j){ 27878569bb4SMatthew G. Knepley ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr); 279a126751eSBarry Smith ierr = PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr); 28078569bb4SMatthew G. Knepley ierr = PetscStrcasecmp(ext_name, var_name, &flg); 28178569bb4SMatthew G. Knepley if (flg) { 28278569bb4SMatthew G. Knepley *varIndex = i+1; 28378569bb4SMatthew G. Knepley PetscFunctionReturn(0); 28478569bb4SMatthew G. Knepley } 28578569bb4SMatthew G. Knepley } 28678569bb4SMatthew G. Knepley } 28778569bb4SMatthew G. Knepley SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unable to locate variable %s in ExodusII file.", name); 28878569bb4SMatthew G. Knepley PetscFunctionReturn(-1); 28978569bb4SMatthew G. Knepley } 29078569bb4SMatthew G. Knepley 29178569bb4SMatthew G. Knepley /* 29278569bb4SMatthew G. Knepley DMPlexView_ExodusII_Internal - Write a DM to disk in exodus format 29378569bb4SMatthew G. Knepley 29478569bb4SMatthew G. Knepley Collective on dm 29578569bb4SMatthew G. Knepley 29678569bb4SMatthew G. Knepley Input Parameters: 29778569bb4SMatthew G. Knepley + dm - The dm to be written 29878569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 29978569bb4SMatthew G. Knepley - degree - the degree of the interpolation space 30078569bb4SMatthew G. Knepley 30178569bb4SMatthew G. Knepley Notes: 30278569bb4SMatthew G. Knepley Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels) 30378569bb4SMatthew G. Knepley consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted. 30478569bb4SMatthew G. Knepley 30578569bb4SMatthew G. Knepley If the dm has been distributed and exoid points to different files on each MPI rank, only the "local" part of the DM 30678569bb4SMatthew G. Knepley (including "ghost" cells and vertices) will be written. If all exoid points to the same file, the resulting file will 30778569bb4SMatthew G. Knepley probably be corrupted. 30878569bb4SMatthew G. Knepley 30978569bb4SMatthew G. Knepley DMPlex only represents geometry while most post-processing software expect that a mesh also provides information 31078569bb4SMatthew G. Knepley on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2. 31178569bb4SMatthew G. Knepley It should be extended to use PetscFE objects. 31278569bb4SMatthew G. Knepley 31378569bb4SMatthew G. Knepley This function will only handle TRI, TET, QUAD and HEX cells. 31478569bb4SMatthew G. Knepley Level: beginner 31578569bb4SMatthew G. Knepley 31678569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadNodal_PlexEXO(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 31778569bb4SMatthew G. Knepley */ 31878569bb4SMatthew G. Knepley PetscErrorCode DMPlexView_ExodusII_Internal(DM dm, int exoid, PetscInt degree) 31978569bb4SMatthew G. Knepley { 32078569bb4SMatthew G. Knepley enum ElemType {TRI, QUAD, TET, HEX}; 32178569bb4SMatthew G. Knepley MPI_Comm comm; 32278569bb4SMatthew G. Knepley /* Connectivity Variables */ 32378569bb4SMatthew G. Knepley PetscInt cellsNotInConnectivity; 32478569bb4SMatthew G. Knepley /* Cell Sets */ 32578569bb4SMatthew G. Knepley DMLabel csLabel; 32678569bb4SMatthew G. Knepley IS csIS; 32778569bb4SMatthew G. Knepley const PetscInt *csIdx; 32878569bb4SMatthew G. Knepley PetscInt num_cs, cs; 32978569bb4SMatthew G. Knepley enum ElemType *type; 330fe209ef9SBlaise Bourdin PetscBool hasLabel; 33178569bb4SMatthew G. Knepley /* Coordinate Variables */ 33278569bb4SMatthew G. Knepley DM cdm; 33378569bb4SMatthew G. Knepley PetscSection section; 33478569bb4SMatthew G. Knepley Vec coord; 33578569bb4SMatthew G. Knepley PetscInt **nodes; 336eae8a387SMatthew G. Knepley PetscInt depth, d, dim, skipCells = 0; 33778569bb4SMatthew G. Knepley PetscInt pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes; 33878569bb4SMatthew G. Knepley PetscInt num_vs, num_fs; 33978569bb4SMatthew G. Knepley PetscMPIInt rank, size; 34078569bb4SMatthew G. Knepley const char *dmName; 341fe209ef9SBlaise Bourdin PetscInt nodesTriP1[4] = {3,0,0,0}; 342fe209ef9SBlaise Bourdin PetscInt nodesTriP2[4] = {3,3,0,0}; 343fe209ef9SBlaise Bourdin PetscInt nodesQuadP1[4] = {4,0,0,0}; 344fe209ef9SBlaise Bourdin PetscInt nodesQuadP2[4] = {4,4,0,1}; 345fe209ef9SBlaise Bourdin PetscInt nodesTetP1[4] = {4,0,0,0}; 346fe209ef9SBlaise Bourdin PetscInt nodesTetP2[4] = {4,6,0,0}; 347fe209ef9SBlaise Bourdin PetscInt nodesHexP1[4] = {8,0,0,0}; 348fe209ef9SBlaise Bourdin PetscInt nodesHexP2[4] = {8,12,6,1}; 34978569bb4SMatthew G. Knepley PetscErrorCode ierr; 35078569bb4SMatthew G. Knepley 35178569bb4SMatthew G. Knepley PetscFunctionBegin; 35278569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 35378569bb4SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 35478569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 35578569bb4SMatthew G. Knepley /* --- Get DM info --- */ 35678569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr); 35778569bb4SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 35878569bb4SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 35978569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 36078569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 36178569bb4SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 36278569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr); 36378569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 36478569bb4SMatthew G. Knepley numCells = cEnd - cStart; 36578569bb4SMatthew G. Knepley numEdges = eEnd - eStart; 36678569bb4SMatthew G. Knepley numVertices = vEnd - vStart; 36778569bb4SMatthew G. Knepley if (depth == 3) {numFaces = fEnd - fStart;} 36878569bb4SMatthew G. Knepley else {numFaces = 0;} 36978569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr); 37078569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr); 37178569bb4SMatthew G. Knepley ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr); 37278569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 37378569bb4SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 37478569bb4SMatthew G. Knepley if (num_cs > 0) { 37578569bb4SMatthew G. Knepley ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr); 37678569bb4SMatthew G. Knepley ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr); 37778569bb4SMatthew G. Knepley ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr); 37878569bb4SMatthew G. Knepley } 37978569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr); 38078569bb4SMatthew G. Knepley /* Set element type for each block and compute total number of nodes */ 38178569bb4SMatthew G. Knepley ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr); 38278569bb4SMatthew G. Knepley numNodes = numVertices; 38378569bb4SMatthew G. Knepley if (degree == 2) {numNodes += numEdges;} 38478569bb4SMatthew G. Knepley cellsNotInConnectivity = numCells; 38578569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 38678569bb4SMatthew G. Knepley IS stratumIS; 38778569bb4SMatthew G. Knepley const PetscInt *cells; 38878569bb4SMatthew G. Knepley PetscScalar *xyz = NULL; 38978569bb4SMatthew G. Knepley PetscInt csSize, closureSize; 39078569bb4SMatthew G. Knepley 39178569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 39278569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 39378569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 39478569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 39578569bb4SMatthew G. Knepley switch (dim) { 39678569bb4SMatthew G. Knepley case 2: 39778569bb4SMatthew G. Knepley if (closureSize == 3*dim) {type[cs] = TRI;} 39878569bb4SMatthew G. Knepley else if (closureSize == 4*dim) {type[cs] = QUAD;} 39978569bb4SMatthew G. Knepley else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 40078569bb4SMatthew G. Knepley break; 40178569bb4SMatthew G. Knepley case 3: 40278569bb4SMatthew G. Knepley if (closureSize == 4*dim) {type[cs] = TET;} 40378569bb4SMatthew G. Knepley else if (closureSize == 8*dim) {type[cs] = HEX;} 40478569bb4SMatthew G. Knepley else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim); 40578569bb4SMatthew G. Knepley break; 40678569bb4SMatthew G. Knepley default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim); 40778569bb4SMatthew G. Knepley } 40878569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;} 40978569bb4SMatthew G. Knepley if ((degree == 2) && (type[cs] == HEX)) {numNodes += csSize; numNodes += numFaces;} 41078569bb4SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr); 41178569bb4SMatthew G. Knepley /* Set nodes and Element type */ 41278569bb4SMatthew G. Knepley if (type[cs] == TRI) { 41378569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTriP1; 41478569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTriP2; 41578569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 41678569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesQuadP1; 41778569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesQuadP2; 41878569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 41978569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesTetP1; 42078569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesTetP2; 42178569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 42278569bb4SMatthew G. Knepley if (degree == 1) nodes[cs] = nodesHexP1; 42378569bb4SMatthew G. Knepley else if (degree == 2) nodes[cs] = nodesHexP2; 42478569bb4SMatthew G. Knepley } 42578569bb4SMatthew G. Knepley /* Compute the number of cells not in the connectivity table */ 42678569bb4SMatthew G. Knepley cellsNotInConnectivity -= nodes[cs][3]*csSize; 42778569bb4SMatthew G. Knepley 42878569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 42978569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 43078569bb4SMatthew G. Knepley } 4316de834b4SMatthew G. Knepley if (num_cs > 0) {PetscStackCallStandard(ex_put_init,(exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs));} 43278569bb4SMatthew G. Knepley /* --- Connectivity --- */ 43378569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 43478569bb4SMatthew G. Knepley IS stratumIS; 43578569bb4SMatthew G. Knepley const PetscInt *cells; 436fe209ef9SBlaise Bourdin PetscInt *connect, off = 0; 437eae8a387SMatthew G. Knepley PetscInt edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0; 43878569bb4SMatthew G. Knepley PetscInt csSize, c, connectSize, closureSize; 43978569bb4SMatthew G. Knepley char *elem_type = NULL; 44078569bb4SMatthew G. Knepley char elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4"; 44178569bb4SMatthew G. Knepley char elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9"; 44278569bb4SMatthew G. Knepley char elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8"; 44378569bb4SMatthew G. Knepley char elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27"; 44478569bb4SMatthew G. Knepley 44578569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 44678569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 44778569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 44878569bb4SMatthew G. Knepley /* Set Element type */ 44978569bb4SMatthew G. Knepley if (type[cs] == TRI) { 45078569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tri3; 45178569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tri6; 45278569bb4SMatthew G. Knepley } else if (type[cs] == QUAD) { 45378569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_quad4; 45478569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_quad9; 45578569bb4SMatthew G. Knepley } else if (type[cs] == TET) { 45678569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_tet4; 45778569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_tet10; 45878569bb4SMatthew G. Knepley } else if (type[cs] == HEX) { 45978569bb4SMatthew G. Knepley if (degree == 1) elem_type = elem_type_hex8; 46078569bb4SMatthew G. Knepley else if (degree == 2) elem_type = elem_type_hex27; 46178569bb4SMatthew G. Knepley } 46278569bb4SMatthew G. Knepley connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3]; 463fe209ef9SBlaise Bourdin ierr = PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect);CHKERRQ(ierr); 464fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_block,(exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1)); 46578569bb4SMatthew G. Knepley /* Find number of vertices, edges, and faces in the closure */ 46678569bb4SMatthew G. Knepley verticesInClosure = nodes[cs][0]; 46778569bb4SMatthew G. Knepley if (depth > 1) { 46878569bb4SMatthew G. Knepley if (dim == 2) { 46978569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr); 47078569bb4SMatthew G. Knepley } else if (dim == 3) { 47178569bb4SMatthew G. Knepley PetscInt *closure = NULL; 47278569bb4SMatthew G. Knepley 47378569bb4SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr); 47478569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 47578569bb4SMatthew G. Knepley edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure; 47678569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 47778569bb4SMatthew G. Knepley } 47878569bb4SMatthew G. Knepley } 47978569bb4SMatthew G. Knepley /* Get connectivity for each cell */ 48078569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 48178569bb4SMatthew G. Knepley PetscInt *closure = NULL; 48278569bb4SMatthew G. Knepley PetscInt temp, i; 48378569bb4SMatthew G. Knepley 48478569bb4SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 48578569bb4SMatthew G. Knepley for (i = 0; i < connectSize; ++i) { 48678569bb4SMatthew G. Knepley if (i < nodes[cs][0]) {/* Vertices */ 487fe209ef9SBlaise Bourdin connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1; 488fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 48978569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */ 490fe209ef9SBlaise Bourdin connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1; 491fe209ef9SBlaise Bourdin if (nodes[cs][2] == 0) connect[i+off] -= numFaces; 492fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 49378569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]){ /* Cells */ 494fe209ef9SBlaise Bourdin connect[i+off] = closure[0] + 1; 495fe209ef9SBlaise Bourdin connect[i+off] -= skipCells; 49678569bb4SMatthew G. Knepley } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]){ /* Faces */ 497fe209ef9SBlaise Bourdin connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1; 498fe209ef9SBlaise Bourdin connect[i+off] -= cellsNotInConnectivity; 49978569bb4SMatthew G. Knepley } else { 500fe209ef9SBlaise Bourdin connect[i+off] = -1; 50178569bb4SMatthew G. Knepley } 50278569bb4SMatthew G. Knepley } 50378569bb4SMatthew G. Knepley /* Tetrahedra are inverted */ 50478569bb4SMatthew G. Knepley if (type[cs] == TET) { 505fe209ef9SBlaise Bourdin temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp; 50678569bb4SMatthew G. Knepley if (degree == 2) { 507e2c9586dSBlaise Bourdin temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp; 508e2c9586dSBlaise Bourdin temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp; 50978569bb4SMatthew G. Knepley } 51078569bb4SMatthew G. Knepley } 51178569bb4SMatthew G. Knepley /* Hexahedra are inverted */ 51278569bb4SMatthew G. Knepley if (type[cs] == HEX) { 513fe209ef9SBlaise Bourdin temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp; 51478569bb4SMatthew G. Knepley if (degree == 2) { 515fe209ef9SBlaise Bourdin temp = connect[8+off]; connect[8+off] = connect[11+off]; connect[11+off] = temp; 516fe209ef9SBlaise Bourdin temp = connect[9+off]; connect[9+off] = connect[10+off]; connect[10+off] = temp; 517fe209ef9SBlaise Bourdin temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp; 518fe209ef9SBlaise Bourdin temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp; 51978569bb4SMatthew G. Knepley 520fe209ef9SBlaise Bourdin temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp; 521fe209ef9SBlaise Bourdin temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp; 522fe209ef9SBlaise Bourdin temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp; 523fe209ef9SBlaise Bourdin temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp; 52478569bb4SMatthew G. Knepley 525fe209ef9SBlaise Bourdin temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp; 526fe209ef9SBlaise Bourdin temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp; 527fe209ef9SBlaise Bourdin temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp; 52878569bb4SMatthew G. Knepley } 52978569bb4SMatthew G. Knepley } 530fe209ef9SBlaise Bourdin off += connectSize; 53178569bb4SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 53278569bb4SMatthew G. Knepley } 533fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_conn,(exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0)); 53478569bb4SMatthew G. Knepley skipCells += (nodes[cs][3] == 0)*csSize; 53578569bb4SMatthew G. Knepley ierr = PetscFree(connect);CHKERRQ(ierr); 53678569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 53778569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 53878569bb4SMatthew G. Knepley } 53978569bb4SMatthew G. Knepley ierr = PetscFree(type);CHKERRQ(ierr); 54078569bb4SMatthew G. Knepley /* --- Coordinates --- */ 54178569bb4SMatthew G. Knepley ierr = PetscSectionCreate(comm, §ion);CHKERRQ(ierr); 54278569bb4SMatthew G. Knepley ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr); 5431e50132fSMatthew G. Knepley if (num_cs) { 54478569bb4SMatthew G. Knepley for (d = 0; d < depth; ++d) { 54578569bb4SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 54678569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 54778569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, p, nodes[0][d] > 0);CHKERRQ(ierr); 54878569bb4SMatthew G. Knepley } 54978569bb4SMatthew G. Knepley } 5501e50132fSMatthew G. Knepley } 55178569bb4SMatthew G. Knepley for (cs = 0; cs < num_cs; ++cs) { 55278569bb4SMatthew G. Knepley IS stratumIS; 55378569bb4SMatthew G. Knepley const PetscInt *cells; 55478569bb4SMatthew G. Knepley PetscInt csSize, c; 55578569bb4SMatthew G. Knepley 55678569bb4SMatthew G. Knepley ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr); 55778569bb4SMatthew G. Knepley ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr); 55878569bb4SMatthew G. Knepley ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr); 55978569bb4SMatthew G. Knepley for (c = 0; c < csSize; ++c) { 56078569bb4SMatthew G. Knepley ierr = PetscSectionSetDof(section, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr); 56178569bb4SMatthew G. Knepley } 56278569bb4SMatthew G. Knepley ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr); 56378569bb4SMatthew G. Knepley ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 56478569bb4SMatthew G. Knepley } 56578569bb4SMatthew G. Knepley if (num_cs > 0) { 56678569bb4SMatthew G. Knepley ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr); 56778569bb4SMatthew G. Knepley ierr = ISDestroy(&csIS);CHKERRQ(ierr); 56878569bb4SMatthew G. Knepley } 56978569bb4SMatthew G. Knepley ierr = PetscFree(nodes);CHKERRQ(ierr); 57078569bb4SMatthew G. Knepley ierr = PetscSectionSetUp(section);CHKERRQ(ierr); 57178569bb4SMatthew G. Knepley if (numNodes > 0) { 57278569bb4SMatthew G. Knepley const char *coordNames[3] = {"x", "y", "z"}; 57378569bb4SMatthew G. Knepley PetscScalar *coords, *closure; 57478569bb4SMatthew G. Knepley PetscReal *cval; 57578569bb4SMatthew G. Knepley PetscInt hasDof, n = 0; 57678569bb4SMatthew G. Knepley 57778569bb4SMatthew G. Knepley /* There can't be more than 24 values in the closure of a point for the coord section */ 5786de834b4SMatthew G. Knepley ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr); 57978569bb4SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr); 58078569bb4SMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 58178569bb4SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 58278569bb4SMatthew G. Knepley ierr = PetscSectionGetDof(section, p, &hasDof); 58378569bb4SMatthew G. Knepley if (hasDof) { 58478569bb4SMatthew G. Knepley PetscInt closureSize = 24, j; 58578569bb4SMatthew G. Knepley 58678569bb4SMatthew G. Knepley ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr); 58778569bb4SMatthew G. Knepley for (d = 0; d < dim; ++d) { 58878569bb4SMatthew G. Knepley cval[d] = 0.0; 58978569bb4SMatthew G. Knepley for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d]; 59078569bb4SMatthew G. Knepley coords[d*numNodes+n] = cval[d] * dim / closureSize; 59178569bb4SMatthew G. Knepley } 59278569bb4SMatthew G. Knepley ++n; 59378569bb4SMatthew G. Knepley } 59478569bb4SMatthew G. Knepley } 5956de834b4SMatthew G. Knepley PetscStackCallStandard(ex_put_coord,(exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes])); 59678569bb4SMatthew G. Knepley ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr); 5976de834b4SMatthew G. Knepley PetscStackCallStandard(ex_put_coord_names,(exoid, (char **) coordNames)); 59878569bb4SMatthew G. Knepley } 59978569bb4SMatthew G. Knepley ierr = PetscSectionDestroy(§ion);CHKERRQ(ierr); 600fe209ef9SBlaise Bourdin /* --- Node Sets/Vertex Sets --- */ 601fe209ef9SBlaise Bourdin DMHasLabel(dm, "Vertex Sets", &hasLabel); 602fe209ef9SBlaise Bourdin if (hasLabel) { 603fe209ef9SBlaise Bourdin PetscInt i, vs, vsSize; 604fe209ef9SBlaise Bourdin const PetscInt *vsIdx, *vertices; 605fe209ef9SBlaise Bourdin PetscInt *nodeList; 606fe209ef9SBlaise Bourdin IS vsIS, stratumIS; 607fe209ef9SBlaise Bourdin DMLabel vsLabel; 608fe209ef9SBlaise Bourdin ierr = DMGetLabel(dm, "Vertex Sets", &vsLabel);CHKERRQ(ierr); 609fe209ef9SBlaise Bourdin ierr = DMLabelGetValueIS(vsLabel, &vsIS);CHKERRQ(ierr); 610fe209ef9SBlaise Bourdin ierr = ISGetIndices(vsIS, &vsIdx);CHKERRQ(ierr); 611fe209ef9SBlaise Bourdin for (vs=0; vs<num_vs; ++vs) { 612fe209ef9SBlaise Bourdin ierr = DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS);CHKERRQ(ierr); 613fe209ef9SBlaise Bourdin ierr = ISGetIndices(stratumIS, &vertices);CHKERRQ(ierr); 614fe209ef9SBlaise Bourdin ierr = ISGetSize(stratumIS, &vsSize);CHKERRQ(ierr); 615fe209ef9SBlaise Bourdin ierr = PetscMalloc1(vsSize, &nodeList); 616fe209ef9SBlaise Bourdin for (i=0; i<vsSize; ++i) { 617fe209ef9SBlaise Bourdin nodeList[i] = vertices[i] - skipCells + 1; 618fe209ef9SBlaise Bourdin } 619fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_set_param,(exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0)); 620fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_set,(exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL)); 621fe209ef9SBlaise Bourdin ierr = ISRestoreIndices(stratumIS, &vertices);CHKERRQ(ierr); 622fe209ef9SBlaise Bourdin ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 623fe209ef9SBlaise Bourdin ierr = PetscFree(nodeList);CHKERRQ(ierr); 624fe209ef9SBlaise Bourdin } 625fe209ef9SBlaise Bourdin ierr = ISRestoreIndices(vsIS, &vsIdx);CHKERRQ(ierr); 626fe209ef9SBlaise Bourdin ierr = ISDestroy(&vsIS);CHKERRQ(ierr); 627fe209ef9SBlaise Bourdin } 628fe209ef9SBlaise Bourdin /* --- Side Sets/Face Sets --- */ 629fe209ef9SBlaise Bourdin ierr = DMHasLabel(dm, "Face Sets", &hasLabel);CHKERRQ(ierr); 630fe209ef9SBlaise Bourdin if (hasLabel) { 631fe209ef9SBlaise Bourdin PetscInt i, j, fs, fsSize; 632fe209ef9SBlaise Bourdin const PetscInt *fsIdx, *faces; 633fe209ef9SBlaise Bourdin IS fsIS, stratumIS; 634fe209ef9SBlaise Bourdin DMLabel fsLabel; 635fe209ef9SBlaise Bourdin PetscInt numPoints, *points; 636fe209ef9SBlaise Bourdin PetscInt elem_list_size = 0; 637fe209ef9SBlaise Bourdin PetscInt *elem_list, *elem_ind, *side_list; 638fe209ef9SBlaise Bourdin 639fe209ef9SBlaise Bourdin ierr = DMGetLabel(dm, "Face Sets", &fsLabel);CHKERRQ(ierr); 640fe209ef9SBlaise Bourdin /* Compute size of Node List and Element List */ 641fe209ef9SBlaise Bourdin ierr = DMLabelGetValueIS(fsLabel, &fsIS);CHKERRQ(ierr); 642fe209ef9SBlaise Bourdin ierr = ISGetIndices(fsIS, &fsIdx);CHKERRQ(ierr); 643fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 644feb237baSPierre Jolivet ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr); 645fe209ef9SBlaise Bourdin ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr); 646fe209ef9SBlaise Bourdin elem_list_size += fsSize; 647fe209ef9SBlaise Bourdin ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 648fe209ef9SBlaise Bourdin } 649fe209ef9SBlaise Bourdin ierr = PetscMalloc3(num_fs, &elem_ind, 650fe209ef9SBlaise Bourdin elem_list_size, &elem_list, 651fe209ef9SBlaise Bourdin elem_list_size, &side_list);CHKERRQ(ierr); 652fe209ef9SBlaise Bourdin elem_ind[0] = 0; 653fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 654fe209ef9SBlaise Bourdin ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr); 655fe209ef9SBlaise Bourdin ierr = ISGetIndices(stratumIS, &faces);CHKERRQ(ierr); 656fe209ef9SBlaise Bourdin ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr); 657fe209ef9SBlaise Bourdin /* Set Parameters */ 658fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_set_param,(exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0)); 659fe209ef9SBlaise Bourdin /* Indices */ 660fe209ef9SBlaise Bourdin if (fs<num_fs-1) { 661fe209ef9SBlaise Bourdin elem_ind[fs+1] = elem_ind[fs] + fsSize; 662fe209ef9SBlaise Bourdin } 663fe209ef9SBlaise Bourdin 664fe209ef9SBlaise Bourdin for (i=0; i<fsSize; ++i) { 665fe209ef9SBlaise Bourdin /* Element List */ 666fe209ef9SBlaise Bourdin points = NULL; 667fe209ef9SBlaise Bourdin ierr = DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr); 668fe209ef9SBlaise Bourdin elem_list[elem_ind[fs] + i] = points[2] +1; 669fe209ef9SBlaise Bourdin ierr = DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr); 670fe209ef9SBlaise Bourdin 671fe209ef9SBlaise Bourdin /* Side List */ 672fe209ef9SBlaise Bourdin points = NULL; 673fe209ef9SBlaise Bourdin ierr = DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 674fe209ef9SBlaise Bourdin for (j=1; j<numPoints; ++j) { 675fe209ef9SBlaise Bourdin if (points[j*2]==faces[i]) {break;} 676fe209ef9SBlaise Bourdin } 677fe209ef9SBlaise Bourdin /* Convert HEX sides */ 678fe209ef9SBlaise Bourdin if (numPoints == 27) { 679fe209ef9SBlaise Bourdin if (j == 1) {j = 5;} 680fe209ef9SBlaise Bourdin else if (j == 2) {j = 6;} 681fe209ef9SBlaise Bourdin else if (j == 3) {j = 1;} 682fe209ef9SBlaise Bourdin else if (j == 4) {j = 3;} 683fe209ef9SBlaise Bourdin else if (j == 5) {j = 2;} 684fe209ef9SBlaise Bourdin else if (j == 6) {j = 4;} 685fe209ef9SBlaise Bourdin } 686fe209ef9SBlaise Bourdin /* Convert TET sides */ 687fe209ef9SBlaise Bourdin if (numPoints == 15) { 688fe209ef9SBlaise Bourdin --j; 689fe209ef9SBlaise Bourdin if (j == 0) {j = 4;} 690fe209ef9SBlaise Bourdin } 691fe209ef9SBlaise Bourdin side_list[elem_ind[fs] + i] = j; 692fe209ef9SBlaise Bourdin ierr = DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr); 693fe209ef9SBlaise Bourdin 694fe209ef9SBlaise Bourdin } 695fe209ef9SBlaise Bourdin ierr = ISRestoreIndices(stratumIS, &faces);CHKERRQ(ierr); 696fe209ef9SBlaise Bourdin ierr = ISDestroy(&stratumIS);CHKERRQ(ierr); 697fe209ef9SBlaise Bourdin } 698fe209ef9SBlaise Bourdin ierr = ISRestoreIndices(fsIS, &fsIdx);CHKERRQ(ierr); 699fe209ef9SBlaise Bourdin ierr = ISDestroy(&fsIS);CHKERRQ(ierr); 700fe209ef9SBlaise Bourdin 701fe209ef9SBlaise Bourdin /* Put side sets */ 702fe209ef9SBlaise Bourdin for (fs=0; fs<num_fs; ++fs) { 703fe209ef9SBlaise Bourdin PetscStackCallStandard(ex_put_set,(exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]])); 704fe209ef9SBlaise Bourdin } 705fe209ef9SBlaise Bourdin ierr = PetscFree3(elem_ind,elem_list,side_list);CHKERRQ(ierr); 706fe209ef9SBlaise Bourdin } 70778569bb4SMatthew G. Knepley PetscFunctionReturn(0); 70878569bb4SMatthew G. Knepley } 70978569bb4SMatthew G. Knepley 71078569bb4SMatthew G. Knepley /* 71178569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file 71278569bb4SMatthew G. Knepley 71378569bb4SMatthew G. Knepley Collective on v 71478569bb4SMatthew G. Knepley 71578569bb4SMatthew G. Knepley Input Parameters: 71678569bb4SMatthew G. Knepley + v - The vector to be written 71778569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 71878569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 71978569bb4SMatthew G. Knepley 72078569bb4SMatthew G. Knepley Notes: 72178569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 72278569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 72378569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 72478569bb4SMatthew G. Knepley amongst all nodal variables. 72578569bb4SMatthew G. Knepley 72678569bb4SMatthew G. Knepley Level: beginner 72778569bb4SMatthew G. Knepley 72878569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO() 72978569bb4SMatthew G. Knepley @*/ 73078569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 73178569bb4SMatthew G. Knepley { 73278569bb4SMatthew G. Knepley MPI_Comm comm; 73378569bb4SMatthew G. Knepley PetscMPIInt size; 73478569bb4SMatthew G. Knepley DM dm; 73578569bb4SMatthew G. Knepley Vec vNatural, vComp; 73622a7544eSVaclav Hapla const PetscScalar *varray; 73778569bb4SMatthew G. Knepley const char *vecname; 73878569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 73978569bb4SMatthew G. Knepley PetscBool useNatural; 74078569bb4SMatthew G. Knepley int offset; 74178569bb4SMatthew G. Knepley PetscErrorCode ierr; 74278569bb4SMatthew G. Knepley 74378569bb4SMatthew G. Knepley PetscFunctionBegin; 74478569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 74578569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 74678569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 74778569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 74878569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 74978569bb4SMatthew G. Knepley if (useNatural) { 75078569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 75178569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 75278569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 75378569bb4SMatthew G. Knepley } else { 75478569bb4SMatthew G. Knepley vNatural = v; 75578569bb4SMatthew G. Knepley } 75678569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 75778569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 75878569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 75978569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file\n", vecname); 76078569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 76178569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 76278569bb4SMatthew G. Knepley We assume that they are stored sequentially */ 76378569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 76478569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 76578569bb4SMatthew G. Knepley if (bs == 1) { 76678569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 7676de834b4SMatthew G. Knepley PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 76878569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 76978569bb4SMatthew G. Knepley } else { 77078569bb4SMatthew G. Knepley IS compIS; 77178569bb4SMatthew G. Knepley PetscInt c; 77278569bb4SMatthew G. Knepley 77378569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 77478569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 77578569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 77678569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 77778569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 7786de834b4SMatthew G. Knepley PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 77978569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 78078569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 78178569bb4SMatthew G. Knepley } 78278569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 78378569bb4SMatthew G. Knepley } 78478569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);} 78578569bb4SMatthew G. Knepley PetscFunctionReturn(0); 78678569bb4SMatthew G. Knepley } 78778569bb4SMatthew G. Knepley 78878569bb4SMatthew G. Knepley /* 78978569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file 79078569bb4SMatthew G. Knepley 79178569bb4SMatthew G. Knepley Collective on v 79278569bb4SMatthew G. Knepley 79378569bb4SMatthew G. Knepley Input Parameters: 79478569bb4SMatthew G. Knepley + v - The vector to be written 79578569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 79678569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 79778569bb4SMatthew G. Knepley 79878569bb4SMatthew G. Knepley Notes: 79978569bb4SMatthew G. Knepley The exodus result nodal variable index is obtained by comparing the Vec name and the 80078569bb4SMatthew G. Knepley names of nodal variables declared in the exodus file. For instance for a Vec named "V" 80178569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 80278569bb4SMatthew G. Knepley amongst all nodal variables. 80378569bb4SMatthew G. Knepley 80478569bb4SMatthew G. Knepley Level: beginner 80578569bb4SMatthew G. Knepley 80678569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO() 80778569bb4SMatthew G. Knepley */ 80878569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step) 80978569bb4SMatthew G. Knepley { 81078569bb4SMatthew G. Knepley MPI_Comm comm; 81178569bb4SMatthew G. Knepley PetscMPIInt size; 81278569bb4SMatthew G. Knepley DM dm; 81378569bb4SMatthew G. Knepley Vec vNatural, vComp; 81422a7544eSVaclav Hapla PetscScalar *varray; 81578569bb4SMatthew G. Knepley const char *vecname; 81678569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 81778569bb4SMatthew G. Knepley PetscBool useNatural; 81878569bb4SMatthew G. Knepley int offset; 81978569bb4SMatthew G. Knepley PetscErrorCode ierr; 82078569bb4SMatthew G. Knepley 82178569bb4SMatthew G. Knepley PetscFunctionBegin; 82278569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 82378569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 82478569bb4SMatthew G. Knepley ierr = VecGetDM(v,&dm);CHKERRQ(ierr); 82578569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 82678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 82778569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 82878569bb4SMatthew G. Knepley else {vNatural = v;} 82978569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 83078569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 83178569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_NODAL, vecname, &offset);CHKERRQ(ierr); 83278569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate nodal variable %s in exodus file", vecname); 83378569bb4SMatthew G. Knepley /* Read local chunk from the file */ 83478569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 83578569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 83678569bb4SMatthew G. Knepley if (bs == 1) { 83778569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 8386de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray)); 83978569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 84078569bb4SMatthew G. Knepley } else { 84178569bb4SMatthew G. Knepley IS compIS; 84278569bb4SMatthew G. Knepley PetscInt c; 84378569bb4SMatthew G. Knepley 84478569bb4SMatthew G. Knepley ierr = ISCreateStride(PETSC_COMM_SELF, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr); 84578569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 84678569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 84778569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 84878569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 8496de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray)); 85078569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 85178569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 85278569bb4SMatthew G. Knepley } 85378569bb4SMatthew G. Knepley ierr = ISDestroy(&compIS);CHKERRQ(ierr); 85478569bb4SMatthew G. Knepley } 85578569bb4SMatthew G. Knepley if (useNatural) { 85678569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 85778569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 85878569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 85978569bb4SMatthew G. Knepley } 86078569bb4SMatthew G. Knepley PetscFunctionReturn(0); 86178569bb4SMatthew G. Knepley } 86278569bb4SMatthew G. Knepley 86378569bb4SMatthew G. Knepley /* 86478569bb4SMatthew G. Knepley VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file 86578569bb4SMatthew G. Knepley 86678569bb4SMatthew G. Knepley Collective on v 86778569bb4SMatthew G. Knepley 86878569bb4SMatthew G. Knepley Input Parameters: 86978569bb4SMatthew G. Knepley + v - The vector to be written 87078569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 87178569bb4SMatthew G. Knepley - step - the time step to write at (exodus steps are numbered starting from 1) 87278569bb4SMatthew G. Knepley 87378569bb4SMatthew G. Knepley Notes: 87478569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 87578569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 87678569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 87778569bb4SMatthew G. Knepley amongst all zonal variables. 87878569bb4SMatthew G. Knepley 87978569bb4SMatthew G. Knepley Level: beginner 88078569bb4SMatthew G. Knepley 88178569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal() 88278569bb4SMatthew G. Knepley */ 88378569bb4SMatthew G. Knepley PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 88478569bb4SMatthew G. Knepley { 88578569bb4SMatthew G. Knepley MPI_Comm comm; 88678569bb4SMatthew G. Knepley PetscMPIInt size; 88778569bb4SMatthew G. Knepley DM dm; 88878569bb4SMatthew G. Knepley Vec vNatural, vComp; 88922a7544eSVaclav Hapla const PetscScalar *varray; 89078569bb4SMatthew G. Knepley const char *vecname; 89178569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 89278569bb4SMatthew G. Knepley PetscBool useNatural; 89378569bb4SMatthew G. Knepley int offset; 89478569bb4SMatthew G. Knepley IS compIS; 89578569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 89678569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 89778569bb4SMatthew G. Knepley PetscErrorCode ierr; 89878569bb4SMatthew G. Knepley 89978569bb4SMatthew G. Knepley PetscFunctionBegin; 90078569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr); 90178569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 90278569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 90378569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 90478569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 90578569bb4SMatthew G. Knepley if (useNatural) { 90678569bb4SMatthew G. Knepley ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr); 90778569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr); 90878569bb4SMatthew G. Knepley ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr); 90978569bb4SMatthew G. Knepley } else { 91078569bb4SMatthew G. Knepley vNatural = v; 91178569bb4SMatthew G. Knepley } 91278569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 91378569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 91478569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 91578569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 91678569bb4SMatthew G. Knepley /* Write local chunk of the result in the exodus file 91778569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 91878569bb4SMatthew G. Knepley We assume that they are stored sequentially 91978569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 92078569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 92178569bb4SMatthew G. Knepley to figure out what to save where. */ 92278569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 92378569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 9246de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 92578569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 92678569bb4SMatthew G. Knepley ex_block block; 92778569bb4SMatthew G. Knepley 92878569bb4SMatthew G. Knepley block.id = csID[set]; 92978569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 9306de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 93178569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 93278569bb4SMatthew G. Knepley } 93378569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 93478569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 93578569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 93678569bb4SMatthew G. Knepley for (set = 0; set < numCS; set++) { 93778569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 93878569bb4SMatthew G. Knepley 93978569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 94078569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 94178569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 94278569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 94378569bb4SMatthew G. Knepley if (bs == 1) { 94478569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr); 9456de834b4SMatthew G. Knepley PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 94678569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr); 94778569bb4SMatthew G. Knepley } else { 94878569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 94978569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 95078569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 95178569bb4SMatthew G. Knepley ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr); 9526de834b4SMatthew G. Knepley 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)])); 95378569bb4SMatthew G. Knepley ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr); 95478569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 95578569bb4SMatthew G. Knepley } 95678569bb4SMatthew G. Knepley } 95778569bb4SMatthew G. Knepley csxs += csSize[set]; 95878569bb4SMatthew G. Knepley } 95978569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 96078569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 96178569bb4SMatthew G. Knepley if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 96278569bb4SMatthew G. Knepley PetscFunctionReturn(0); 96378569bb4SMatthew G. Knepley } 96478569bb4SMatthew G. Knepley 96578569bb4SMatthew G. Knepley /* 96678569bb4SMatthew G. Knepley VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file 96778569bb4SMatthew G. Knepley 96878569bb4SMatthew G. Knepley Collective on v 96978569bb4SMatthew G. Knepley 97078569bb4SMatthew G. Knepley Input Parameters: 97178569bb4SMatthew G. Knepley + v - The vector to be written 97278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance) 97378569bb4SMatthew G. Knepley - step - the time step to read at (exodus steps are numbered starting from 1) 97478569bb4SMatthew G. Knepley 97578569bb4SMatthew G. Knepley Notes: 97678569bb4SMatthew G. Knepley The exodus result zonal variable index is obtained by comparing the Vec name and the 97778569bb4SMatthew G. Knepley names of zonal variables declared in the exodus file. For instance for a Vec named "V" 97878569bb4SMatthew G. Knepley the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11" 97978569bb4SMatthew G. Knepley amongst all zonal variables. 98078569bb4SMatthew G. Knepley 98178569bb4SMatthew G. Knepley Level: beginner 98278569bb4SMatthew G. Knepley 98378569bb4SMatthew G. Knepley .seealso: EXOGetVarIndex_Private(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal() 98478569bb4SMatthew G. Knepley */ 98578569bb4SMatthew G. Knepley PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step) 98678569bb4SMatthew G. Knepley { 98778569bb4SMatthew G. Knepley MPI_Comm comm; 98878569bb4SMatthew G. Knepley PetscMPIInt size; 98978569bb4SMatthew G. Knepley DM dm; 99078569bb4SMatthew G. Knepley Vec vNatural, vComp; 99122a7544eSVaclav Hapla PetscScalar *varray; 99278569bb4SMatthew G. Knepley const char *vecname; 99378569bb4SMatthew G. Knepley PetscInt xs, xe, bs; 99478569bb4SMatthew G. Knepley PetscBool useNatural; 99578569bb4SMatthew G. Knepley int offset; 99678569bb4SMatthew G. Knepley IS compIS; 99778569bb4SMatthew G. Knepley PetscInt *csSize, *csID; 99878569bb4SMatthew G. Knepley PetscInt numCS, set, csxs = 0; 99978569bb4SMatthew G. Knepley PetscErrorCode ierr; 100078569bb4SMatthew G. Knepley 100178569bb4SMatthew G. Knepley PetscFunctionBegin; 100278569bb4SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr); 100378569bb4SMatthew G. Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 100478569bb4SMatthew G. Knepley ierr = VecGetDM(v, &dm);CHKERRQ(ierr); 100578569bb4SMatthew G. Knepley ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr); 100678569bb4SMatthew G. Knepley useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE; 100778569bb4SMatthew G. Knepley if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);} 100878569bb4SMatthew G. Knepley else {vNatural = v;} 100978569bb4SMatthew G. Knepley /* Get the location of the variable in the exodus file */ 101078569bb4SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr); 101178569bb4SMatthew G. Knepley ierr = EXOGetVarIndex_Private(exoid, EX_ELEM_BLOCK, vecname, &offset);CHKERRQ(ierr); 101278569bb4SMatthew G. Knepley if (!offset) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Cannot locate zonal variable %s in exodus file\n", vecname); 101378569bb4SMatthew G. Knepley /* Read local chunk of the result in the exodus file 101478569bb4SMatthew G. Knepley exodus stores each component of a vector-valued field as a separate variable. 101578569bb4SMatthew G. Knepley We assume that they are stored sequentially 101678569bb4SMatthew G. Knepley Zonal variables are accessed one element block at a time, so we loop through the cell sets, 101778569bb4SMatthew G. Knepley but once the vector has been reordered to natural size, we cannot use the label informations 101878569bb4SMatthew G. Knepley to figure out what to save where. */ 101978569bb4SMatthew G. Knepley numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK); 102078569bb4SMatthew G. Knepley ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr); 10216de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID)); 102278569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 102378569bb4SMatthew G. Knepley ex_block block; 102478569bb4SMatthew G. Knepley 102578569bb4SMatthew G. Knepley block.id = csID[set]; 102678569bb4SMatthew G. Knepley block.type = EX_ELEM_BLOCK; 10276de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_block_param,(exoid, &block)); 102878569bb4SMatthew G. Knepley csSize[set] = block.num_entry; 102978569bb4SMatthew G. Knepley } 103078569bb4SMatthew G. Knepley ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr); 103178569bb4SMatthew G. Knepley ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr); 103278569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);} 103378569bb4SMatthew G. Knepley for (set = 0; set < numCS; ++set) { 103478569bb4SMatthew G. Knepley PetscInt csLocalSize, c; 103578569bb4SMatthew G. Knepley 103678569bb4SMatthew G. Knepley /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1 103778569bb4SMatthew G. Knepley local slice of zonal values: xs/bs,xm/bs-1 103878569bb4SMatthew G. Knepley intersection: max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */ 103978569bb4SMatthew G. Knepley csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs)); 104078569bb4SMatthew G. Knepley if (bs == 1) { 104178569bb4SMatthew G. Knepley ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr); 10426de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)])); 104378569bb4SMatthew G. Knepley ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr); 104478569bb4SMatthew G. Knepley } else { 104578569bb4SMatthew G. Knepley for (c = 0; c < bs; ++c) { 104678569bb4SMatthew G. Knepley ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr); 104778569bb4SMatthew G. Knepley ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 104878569bb4SMatthew G. Knepley ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr); 10496de834b4SMatthew G. Knepley 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)])); 105078569bb4SMatthew G. Knepley ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr); 105178569bb4SMatthew G. Knepley ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr); 105278569bb4SMatthew G. Knepley } 105378569bb4SMatthew G. Knepley } 105478569bb4SMatthew G. Knepley csxs += csSize[set]; 105578569bb4SMatthew G. Knepley } 105678569bb4SMatthew G. Knepley ierr = PetscFree2(csID, csSize);CHKERRQ(ierr); 105778569bb4SMatthew G. Knepley if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);} 105878569bb4SMatthew G. Knepley if (useNatural) { 105978569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr); 106078569bb4SMatthew G. Knepley ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr); 106178569bb4SMatthew G. Knepley ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr); 106278569bb4SMatthew G. Knepley } 106378569bb4SMatthew G. Knepley PetscFunctionReturn(0); 106478569bb4SMatthew G. Knepley } 1065b53c8456SSatish Balay #endif 106678569bb4SMatthew G. Knepley 10671e50132fSMatthew G. Knepley /*@ 10681e50132fSMatthew G. Knepley PetscViewerExodusIIGetId - Get the file id of the ExodusII file 10691e50132fSMatthew G. Knepley 10701e50132fSMatthew G. Knepley Logically Collective on PetscViewer 10711e50132fSMatthew G. Knepley 10721e50132fSMatthew G. Knepley Input Parameter: 10731e50132fSMatthew G. Knepley . viewer - the PetscViewer 10741e50132fSMatthew G. Knepley 10751e50132fSMatthew G. Knepley Output Parameter: 10761e50132fSMatthew G. Knepley - exoid - The ExodusII file id 10771e50132fSMatthew G. Knepley 10781e50132fSMatthew G. Knepley Level: intermediate 10791e50132fSMatthew G. Knepley 10801e50132fSMatthew G. Knepley .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 10811e50132fSMatthew G. Knepley @*/ 10821e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) 10831e50132fSMatthew G. Knepley { 10841e50132fSMatthew G. Knepley PetscErrorCode ierr; 10851e50132fSMatthew G. Knepley 10861e50132fSMatthew G. Knepley PetscFunctionBegin; 10871e50132fSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 10881e50132fSMatthew G. Knepley ierr = PetscTryMethod(viewer, "PetscViewerExodusIIGetId_C",(PetscViewer,int*),(viewer,exoid));CHKERRQ(ierr); 10891e50132fSMatthew G. Knepley PetscFunctionReturn(0); 10901e50132fSMatthew G. Knepley } 10911e50132fSMatthew G. Knepley 109233751fbdSMichael Lange /*@C 1093eaf898f9SPatrick Sanan DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file. 109433751fbdSMichael Lange 1095d083f849SBarry Smith Collective 109633751fbdSMichael Lange 109733751fbdSMichael Lange Input Parameters: 109833751fbdSMichael Lange + comm - The MPI communicator 109933751fbdSMichael Lange . filename - The name of the ExodusII file 110033751fbdSMichael Lange - interpolate - Create faces and edges in the mesh 110133751fbdSMichael Lange 110233751fbdSMichael Lange Output Parameter: 110333751fbdSMichael Lange . dm - The DM object representing the mesh 110433751fbdSMichael Lange 110533751fbdSMichael Lange Level: beginner 110633751fbdSMichael Lange 110733751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus() 110833751fbdSMichael Lange @*/ 110933751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) 111033751fbdSMichael Lange { 111133751fbdSMichael Lange PetscMPIInt rank; 111233751fbdSMichael Lange PetscErrorCode ierr; 111333751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 1114e84f7031SBlaise Bourdin int CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1; 111533751fbdSMichael Lange float version; 111633751fbdSMichael Lange #endif 111733751fbdSMichael Lange 111833751fbdSMichael Lange PetscFunctionBegin; 111933751fbdSMichael Lange PetscValidCharPointer(filename, 2); 112033751fbdSMichael Lange ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 112133751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII) 112233751fbdSMichael Lange if (!rank) { 112333751fbdSMichael Lange exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version); 112433751fbdSMichael Lange if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename); 112533751fbdSMichael Lange } 112633751fbdSMichael Lange ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr); 11276de834b4SMatthew G. Knepley if (!rank) {PetscStackCallStandard(ex_close,(exoid));} 112833751fbdSMichael Lange #else 112933751fbdSMichael Lange SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 113033751fbdSMichael Lange #endif 113133751fbdSMichael Lange PetscFunctionReturn(0); 113233751fbdSMichael Lange } 113333751fbdSMichael Lange 1134552f7358SJed Brown /*@ 113533751fbdSMichael Lange DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID. 1136552f7358SJed Brown 1137d083f849SBarry Smith Collective 1138552f7358SJed Brown 1139552f7358SJed Brown Input Parameters: 1140552f7358SJed Brown + comm - The MPI communicator 1141552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open 1142552f7358SJed Brown - interpolate - Create faces and edges in the mesh 1143552f7358SJed Brown 1144552f7358SJed Brown Output Parameter: 1145552f7358SJed Brown . dm - The DM object representing the mesh 1146552f7358SJed Brown 1147552f7358SJed Brown Level: beginner 1148552f7358SJed Brown 1149e84e7769SJed Brown .seealso: DMPLEX, DMCreate() 1150552f7358SJed Brown @*/ 1151552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) 1152552f7358SJed Brown { 1153552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 1154552f7358SJed Brown PetscMPIInt num_proc, rank; 1155552f7358SJed Brown PetscSection coordSection; 1156552f7358SJed Brown Vec coordinates; 1157552f7358SJed Brown PetscScalar *coords; 1158552f7358SJed Brown PetscInt coordSize, v; 1159552f7358SJed Brown PetscErrorCode ierr; 1160552f7358SJed Brown /* Read from ex_get_init() */ 1161552f7358SJed Brown char title[PETSC_MAX_PATH_LEN+1]; 1162e8f6893fSMatthew G. Knepley int dim = 0, numVertices = 0, numCells = 0, numHybridCells = 0; 1163552f7358SJed Brown int num_cs = 0, num_vs = 0, num_fs = 0; 1164552f7358SJed Brown #endif 1165552f7358SJed Brown 1166552f7358SJed Brown PetscFunctionBegin; 1167552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII) 1168552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1169552f7358SJed Brown ierr = MPI_Comm_size(comm, &num_proc);CHKERRQ(ierr); 1170552f7358SJed Brown ierr = DMCreate(comm, dm);CHKERRQ(ierr); 1171552f7358SJed Brown ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr); 1172552f7358SJed Brown /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */ 1173552f7358SJed Brown if (!rank) { 1174580bdb30SBarry Smith ierr = PetscMemzero(title,PETSC_MAX_PATH_LEN+1);CHKERRQ(ierr); 11756de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_init,(exoid, title, &dim, &numVertices, &numCells, &num_cs, &num_vs, &num_fs)); 1176552f7358SJed Brown if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n"); 1177552f7358SJed Brown } 1178552f7358SJed Brown ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRQ(ierr); 1179552f7358SJed Brown ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRQ(ierr); 1180552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr); 1181c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dm, dim);CHKERRQ(ierr); 1182552f7358SJed Brown ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr); 1183*412e9a14SMatthew G. Knepley /* We do not want this label automatically computed, instead we compute it here */ 1184*412e9a14SMatthew G. Knepley ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr); 1185552f7358SJed Brown 1186552f7358SJed Brown /* Read cell sets information */ 1187552f7358SJed Brown if (!rank) { 1188552f7358SJed Brown PetscInt *cone; 1189e8f6893fSMatthew G. Knepley int c, cs, ncs, c_loc, v, v_loc; 1190552f7358SJed Brown /* Read from ex_get_elem_blk_ids() */ 1191e8f6893fSMatthew G. Knepley int *cs_id, *cs_order; 1192552f7358SJed Brown /* Read from ex_get_elem_block() */ 1193552f7358SJed Brown char buffer[PETSC_MAX_PATH_LEN+1]; 1194e8f6893fSMatthew G. Knepley int num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr; 1195552f7358SJed Brown /* Read from ex_get_elem_conn() */ 1196552f7358SJed Brown int *cs_connect; 1197552f7358SJed Brown 1198552f7358SJed Brown /* Get cell sets IDs */ 1199e8f6893fSMatthew G. Knepley ierr = PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order);CHKERRQ(ierr); 12006de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id)); 1201552f7358SJed Brown /* Read the cell set connectivity table and build mesh topology 1202552f7358SJed Brown EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */ 1203e8f6893fSMatthew G. Knepley /* Check for a hybrid mesh */ 1204e8f6893fSMatthew G. Knepley for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) { 1205e8f6893fSMatthew G. Knepley PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr)); 1206e8f6893fSMatthew G. Knepley switch (dim) { 1207e8f6893fSMatthew G. Knepley case 3: 1208e8f6893fSMatthew G. Knepley switch (num_vertex_per_cell) { 1209e8f6893fSMatthew G. Knepley case 6: 1210e8f6893fSMatthew G. Knepley cs_order[cs] = cs; 1211e8f6893fSMatthew G. Knepley numHybridCells += num_cell_in_set; 1212e8f6893fSMatthew G. Knepley ++num_hybrid; 1213e8f6893fSMatthew G. Knepley break; 1214e8f6893fSMatthew G. Knepley default: 1215e8f6893fSMatthew G. Knepley for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1]; 1216e8f6893fSMatthew G. Knepley cs_order[cs-num_hybrid] = cs; 1217e8f6893fSMatthew G. Knepley } 1218e8f6893fSMatthew G. Knepley break; 1219e8f6893fSMatthew G. Knepley default: 1220e8f6893fSMatthew G. Knepley for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1]; 1221e8f6893fSMatthew G. Knepley cs_order[cs-num_hybrid] = cs; 1222e8f6893fSMatthew G. Knepley } 1223e8f6893fSMatthew G. Knepley } 1224552f7358SJed Brown /* First set sizes */ 1225e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1226e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 12276de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr)); 1228552f7358SJed Brown for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1229552f7358SJed Brown ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr); 1230*412e9a14SMatthew G. Knepley if (c >= numCells-numHybridCells) { 1231*412e9a14SMatthew G. Knepley ierr = DMPlexSetCellType(*dm, c, DM_POLYTOPE_TRI_PRISM);CHKERRQ(ierr); 1232*412e9a14SMatthew G. Knepley } else { 1233*412e9a14SMatthew G. Knepley DMPolytopeType ct; 1234*412e9a14SMatthew G. Knepley 1235*412e9a14SMatthew G. Knepley ierr = DMPlexComputeCellType_Internal(*dm, c, 1, &ct);CHKERRQ(ierr); 1236*412e9a14SMatthew G. Knepley ierr = DMPlexSetCellType(*dm, c, ct);CHKERRQ(ierr); 1237552f7358SJed Brown } 1238552f7358SJed Brown } 1239*412e9a14SMatthew G. Knepley } 1240*412e9a14SMatthew G. Knepley for (v = numCells; v < numCells+numVertices; ++v) {ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr);} 1241552f7358SJed Brown ierr = DMSetUp(*dm);CHKERRQ(ierr); 1242e8f6893fSMatthew G. Knepley for (ncs = 0, c = 0; ncs < num_cs; ++ncs) { 1243e8f6893fSMatthew G. Knepley const PetscInt cs = cs_order[ncs]; 12446de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr)); 1245dcca6d9dSJed Brown ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr); 12466de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL)); 1247eaf898f9SPatrick Sanan /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */ 1248552f7358SJed Brown for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) { 1249552f7358SJed Brown for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) { 1250552f7358SJed Brown cone[v_loc] = cs_connect[v]+numCells-1; 1251552f7358SJed Brown } 1252731dcdf4SMatthew G. Knepley if (dim == 3) { 12532e1b13c2SMatthew G. Knepley /* Tetrahedra are inverted */ 12542e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 4) { 12552e1b13c2SMatthew G. Knepley PetscInt tmp = cone[0]; 12562e1b13c2SMatthew G. Knepley cone[0] = cone[1]; 12572e1b13c2SMatthew G. Knepley cone[1] = tmp; 12582e1b13c2SMatthew G. Knepley } 12592e1b13c2SMatthew G. Knepley /* Hexahedra are inverted */ 12602e1b13c2SMatthew G. Knepley if (num_vertex_per_cell == 8) { 12612e1b13c2SMatthew G. Knepley PetscInt tmp = cone[1]; 12622e1b13c2SMatthew G. Knepley cone[1] = cone[3]; 12632e1b13c2SMatthew G. Knepley cone[3] = tmp; 12642e1b13c2SMatthew G. Knepley } 1265*412e9a14SMatthew G. Knepley /* Triangular prisms are inverted */ 1266*412e9a14SMatthew G. Knepley if (num_vertex_per_cell == 6) { 1267*412e9a14SMatthew G. Knepley PetscInt tmp = cone[1]; 1268*412e9a14SMatthew G. Knepley cone[1] = cone[2]; 1269*412e9a14SMatthew G. Knepley cone[2] = tmp; 1270*412e9a14SMatthew G. Knepley } 1271731dcdf4SMatthew G. Knepley } 1272552f7358SJed Brown ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr); 1273c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr); 1274552f7358SJed Brown } 1275552f7358SJed Brown ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr); 1276552f7358SJed Brown } 1277e8f6893fSMatthew G. Knepley ierr = PetscFree2(cs_id, cs_order);CHKERRQ(ierr); 1278552f7358SJed Brown } 1279552f7358SJed Brown ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr); 1280552f7358SJed Brown ierr = DMPlexStratify(*dm);CHKERRQ(ierr); 1281552f7358SJed Brown if (interpolate) { 12825fd9971aSMatthew G. Knepley DM idm; 1283552f7358SJed Brown 12849f074e33SMatthew G Knepley ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr); 1285552f7358SJed Brown ierr = DMDestroy(dm);CHKERRQ(ierr); 1286552f7358SJed Brown *dm = idm; 1287552f7358SJed Brown } 1288552f7358SJed Brown 1289552f7358SJed Brown /* Create vertex set label */ 1290552f7358SJed Brown if (!rank && (num_vs > 0)) { 1291552f7358SJed Brown int vs, v; 1292552f7358SJed Brown /* Read from ex_get_node_set_ids() */ 1293552f7358SJed Brown int *vs_id; 1294552f7358SJed Brown /* Read from ex_get_node_set_param() */ 1295f958083aSBlaise Bourdin int num_vertex_in_set; 1296552f7358SJed Brown /* Read from ex_get_node_set() */ 1297552f7358SJed Brown int *vs_vertex_list; 1298552f7358SJed Brown 1299552f7358SJed Brown /* Get vertex set ids */ 1300785e854fSJed Brown ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr); 13016de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id)); 1302552f7358SJed Brown for (vs = 0; vs < num_vs; ++vs) { 13036de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL)); 1304785e854fSJed Brown ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr); 13056de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL)); 1306552f7358SJed Brown for (v = 0; v < num_vertex_in_set; ++v) { 1307c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr); 1308552f7358SJed Brown } 1309552f7358SJed Brown ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr); 1310552f7358SJed Brown } 1311552f7358SJed Brown ierr = PetscFree(vs_id);CHKERRQ(ierr); 1312552f7358SJed Brown } 1313552f7358SJed Brown /* Read coordinates */ 131469d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr); 1315552f7358SJed Brown ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr); 1316552f7358SJed Brown ierr = PetscSectionSetFieldComponents(coordSection, 0, dim);CHKERRQ(ierr); 1317552f7358SJed Brown ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr); 1318552f7358SJed Brown for (v = numCells; v < numCells+numVertices; ++v) { 1319552f7358SJed Brown ierr = PetscSectionSetDof(coordSection, v, dim);CHKERRQ(ierr); 1320552f7358SJed Brown ierr = PetscSectionSetFieldDof(coordSection, v, 0, dim);CHKERRQ(ierr); 1321552f7358SJed Brown } 1322552f7358SJed Brown ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr); 1323552f7358SJed Brown ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr); 13248b9ced59SLisandro Dalcin ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr); 1325552f7358SJed Brown ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr); 1326552f7358SJed Brown ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr); 13274e90ef8eSMatthew G. Knepley ierr = VecSetBlockSize(coordinates, dim);CHKERRQ(ierr); 13282eb5907fSJed Brown ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr); 1329552f7358SJed Brown ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 13300aba08f6SKarl Rupp if (!rank) { 1331e84f7031SBlaise Bourdin PetscReal *x, *y, *z; 1332552f7358SJed Brown 1333dcca6d9dSJed Brown ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr); 13346de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_coord,(exoid, x, y, z)); 13350d644c17SKarl Rupp if (dim > 0) { 13360d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+0] = x[v]; 13370d644c17SKarl Rupp } 13380d644c17SKarl Rupp if (dim > 1) { 13390d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+1] = y[v]; 13400d644c17SKarl Rupp } 13410d644c17SKarl Rupp if (dim > 2) { 13420d644c17SKarl Rupp for (v = 0; v < numVertices; ++v) coords[v*dim+2] = z[v]; 13430d644c17SKarl Rupp } 1344552f7358SJed Brown ierr = PetscFree3(x,y,z);CHKERRQ(ierr); 1345552f7358SJed Brown } 1346552f7358SJed Brown ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 1347552f7358SJed Brown ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr); 1348552f7358SJed Brown ierr = VecDestroy(&coordinates);CHKERRQ(ierr); 1349552f7358SJed Brown 1350552f7358SJed Brown /* Create side set label */ 1351552f7358SJed Brown if (!rank && interpolate && (num_fs > 0)) { 1352552f7358SJed Brown int fs, f, voff; 1353552f7358SJed Brown /* Read from ex_get_side_set_ids() */ 1354552f7358SJed Brown int *fs_id; 1355552f7358SJed Brown /* Read from ex_get_side_set_param() */ 1356f958083aSBlaise Bourdin int num_side_in_set; 1357552f7358SJed Brown /* Read from ex_get_side_set_node_list() */ 1358552f7358SJed Brown int *fs_vertex_count_list, *fs_vertex_list; 1359ef073283Smichael_afanasiev /* Read side set labels */ 1360ef073283Smichael_afanasiev char fs_name[MAX_STR_LENGTH+1]; 1361552f7358SJed Brown 1362552f7358SJed Brown /* Get side set ids */ 1363785e854fSJed Brown ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr); 13646de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id)); 1365552f7358SJed Brown for (fs = 0; fs < num_fs; ++fs) { 13666de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL)); 1367dcca6d9dSJed Brown ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr); 13686de834b4SMatthew G. Knepley PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list)); 1369ef073283Smichael_afanasiev /* Get the specific name associated with this side set ID. */ 1370ef073283Smichael_afanasiev int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name); 1371552f7358SJed Brown for (f = 0, voff = 0; f < num_side_in_set; ++f) { 13720298fd71SBarry Smith const PetscInt *faces = NULL; 1373552f7358SJed Brown PetscInt faceSize = fs_vertex_count_list[f], numFaces; 1374552f7358SJed Brown PetscInt faceVertices[4], v; 1375552f7358SJed Brown 1376552f7358SJed Brown if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize); 1377552f7358SJed Brown for (v = 0; v < faceSize; ++v, ++voff) { 1378552f7358SJed Brown faceVertices[v] = fs_vertex_list[voff]+numCells-1; 1379552f7358SJed Brown } 1380552f7358SJed Brown ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 1381552f7358SJed Brown if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces); 1382c58f1c22SToby Isaac ierr = DMSetLabelValue(*dm, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr); 1383ef073283Smichael_afanasiev /* Only add the label if one has been detected for this side set. */ 1384ef073283Smichael_afanasiev if (!fs_name_err) { 1385ef073283Smichael_afanasiev ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr); 1386ef073283Smichael_afanasiev } 1387552f7358SJed Brown ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr); 1388552f7358SJed Brown } 1389552f7358SJed Brown ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr); 1390552f7358SJed Brown } 1391552f7358SJed Brown ierr = PetscFree(fs_id);CHKERRQ(ierr); 1392552f7358SJed Brown } 1393552f7358SJed Brown #else 1394552f7358SJed Brown SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii"); 1395552f7358SJed Brown #endif 1396552f7358SJed Brown PetscFunctionReturn(0); 1397552f7358SJed Brown } 1398