xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown 
4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII)
5552f7358SJed Brown #include <netcdf.h>
6552f7358SJed Brown #include <exodusII.h>
7552f7358SJed Brown #endif
8552f7358SJed Brown 
91e50132fSMatthew G. Knepley #include <petsc/private/viewerimpl.h>
101e50132fSMatthew G. Knepley #include <petsc/private/viewerexodusiiimpl.h>
11e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
1278569bb4SMatthew G. Knepley /*
131e50132fSMatthew G. Knepley   PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator.
141e50132fSMatthew G. Knepley 
151e50132fSMatthew G. Knepley   Collective
161e50132fSMatthew G. Knepley 
171e50132fSMatthew G. Knepley   Input Parameter:
181e50132fSMatthew G. Knepley . comm - the MPI communicator to share the ExodusII PetscViewer
191e50132fSMatthew G. Knepley 
201e50132fSMatthew G. Knepley   Level: intermediate
211e50132fSMatthew G. Knepley 
221e50132fSMatthew G. Knepley   Notes:
231e50132fSMatthew G. Knepley     misses Fortran bindings
241e50132fSMatthew G. Knepley 
251e50132fSMatthew G. Knepley   Notes:
261e50132fSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return
271e50132fSMatthew G. Knepley   an error code.  The GLVIS PetscViewer is usually used in the form
281e50132fSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
291e50132fSMatthew G. Knepley 
30db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()`
311e50132fSMatthew G. Knepley */
329371c9d4SSatish Balay PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm) {
331e50132fSMatthew G. Knepley   PetscViewer    viewer;
341e50132fSMatthew G. Knepley   PetscErrorCode ierr;
351e50132fSMatthew G. Knepley 
361e50132fSMatthew G. Knepley   PetscFunctionBegin;
371e50132fSMatthew G. Knepley   ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer);
389371c9d4SSatish Balay   if (ierr) {
399371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
4011cc89d2SBarry Smith     PetscFunctionReturn(NULL);
419371c9d4SSatish Balay   }
421e50132fSMatthew G. Knepley   ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
439371c9d4SSatish Balay   if (ierr) {
449371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
4511cc89d2SBarry Smith     PetscFunctionReturn(NULL);
469371c9d4SSatish Balay   }
471e50132fSMatthew G. Knepley   PetscFunctionReturn(viewer);
481e50132fSMatthew G. Knepley }
491e50132fSMatthew G. Knepley 
509371c9d4SSatish Balay static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer) {
516823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data;
521e50132fSMatthew G. Knepley 
531e50132fSMatthew G. Knepley   PetscFunctionBegin;
549566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename:    %s\n", exo->filename));
559566063dSJacob Faibussowitsch   if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid:       %d\n", exo->exoid));
569566063dSJacob Faibussowitsch   if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode:     %d\n", exo->btype));
57197f6770SJed Brown   if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order:  %" PetscInt_FMT "\n", exo->order));
581e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
591e50132fSMatthew G. Knepley }
601e50132fSMatthew G. Knepley 
619371c9d4SSatish Balay static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscViewer v, PetscOptionItems *PetscOptionsObject) {
621e50132fSMatthew G. Knepley   PetscFunctionBegin;
63d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options");
64d0609cedSBarry Smith   PetscOptionsHeadEnd();
651e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
661e50132fSMatthew G. Knepley }
671e50132fSMatthew G. Knepley 
689371c9d4SSatish Balay static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer) {
691e50132fSMatthew G. Knepley   PetscFunctionBegin;
701e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
711e50132fSMatthew G. Knepley }
721e50132fSMatthew G. Knepley 
739371c9d4SSatish Balay static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer) {
741e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
751e50132fSMatthew G. Knepley 
761e50132fSMatthew G. Knepley   PetscFunctionBegin;
7748a46eb9SPierre Jolivet   if (exo->exoid >= 0) PetscCallExternal(ex_close, exo->exoid);
789566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo->filename));
799566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo));
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
832e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetId_C", NULL));
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetOrder_C", NULL));
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerSetOrder_C", NULL));
871e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
881e50132fSMatthew G. Knepley }
891e50132fSMatthew G. Knepley 
909371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[]) {
911e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
921e50132fSMatthew G. Knepley   PetscMPIInt           rank;
931e50132fSMatthew G. Knepley   int                   CPU_word_size, IO_word_size, EXO_mode;
946823f3c5SBlaise Bourdin   MPI_Info              mpi_info = MPI_INFO_NULL;
956823f3c5SBlaise Bourdin   float                 EXO_version;
961e50132fSMatthew G. Knepley 
979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
981e50132fSMatthew G. Knepley   CPU_word_size = sizeof(PetscReal);
991e50132fSMatthew G. Knepley   IO_word_size  = sizeof(PetscReal);
1001e50132fSMatthew G. Knepley 
1011e50132fSMatthew G. Knepley   PetscFunctionBegin;
1026823f3c5SBlaise Bourdin   if (exo->exoid >= 0) {
103792fecdfSBarry Smith     PetscCallExternal(ex_close, exo->exoid);
1046823f3c5SBlaise Bourdin     exo->exoid = -1;
1056823f3c5SBlaise Bourdin   }
1069566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscFree(exo->filename));
1079566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &exo->filename));
1081e50132fSMatthew G. Knepley   switch (exo->btype) {
1099371c9d4SSatish Balay   case FILE_MODE_READ: EXO_mode = EX_READ; break;
1101e50132fSMatthew G. Knepley   case FILE_MODE_APPEND:
1116823f3c5SBlaise Bourdin   case FILE_MODE_UPDATE:
1126823f3c5SBlaise Bourdin   case FILE_MODE_APPEND_UPDATE:
1136823f3c5SBlaise Bourdin     /* Will fail if the file does not already exist */
1146823f3c5SBlaise Bourdin     EXO_mode = EX_WRITE;
1151e50132fSMatthew G. Knepley     break;
1161e50132fSMatthew G. Knepley   case FILE_MODE_WRITE:
1176823f3c5SBlaise Bourdin     /*
1186823f3c5SBlaise Bourdin       exodus only allows writing geometry upon file creation, so we will let DMView create the file.
1196823f3c5SBlaise Bourdin     */
1206823f3c5SBlaise Bourdin     PetscFunctionReturn(0);
1211e50132fSMatthew G. Knepley     break;
1229371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
1231e50132fSMatthew G. Knepley   }
1241e50132fSMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
1251e50132fSMatthew G. Knepley   EXO_mode += EX_ALL_INT64_API;
1261e50132fSMatthew G. Knepley #endif
1276823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(name, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PETSC_COMM_WORLD, mpi_info);
12808401ef6SPierre Jolivet   PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
1291e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1301e50132fSMatthew G. Knepley }
1311e50132fSMatthew G. Knepley 
1329371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name) {
1331e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1341e50132fSMatthew G. Knepley 
1351e50132fSMatthew G. Knepley   PetscFunctionBegin;
1361e50132fSMatthew G. Knepley   *name = exo->filename;
1371e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1381e50132fSMatthew G. Knepley }
1391e50132fSMatthew G. Knepley 
1409371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type) {
1411e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1421e50132fSMatthew G. Knepley 
1431e50132fSMatthew G. Knepley   PetscFunctionBegin;
1441e50132fSMatthew G. Knepley   exo->btype = type;
1451e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1461e50132fSMatthew G. Knepley }
1471e50132fSMatthew G. Knepley 
1489371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type) {
1491e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1501e50132fSMatthew G. Knepley 
1511e50132fSMatthew G. Knepley   PetscFunctionBegin;
1521e50132fSMatthew G. Knepley   *type = exo->btype;
1531e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1541e50132fSMatthew G. Knepley }
1551e50132fSMatthew G. Knepley 
1569371c9d4SSatish Balay static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid) {
1571e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1581e50132fSMatthew G. Knepley 
1591e50132fSMatthew G. Knepley   PetscFunctionBegin;
1601e50132fSMatthew G. Knepley   *exoid = exo->exoid;
1611e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1621e50132fSMatthew G. Knepley }
1631e50132fSMatthew G. Knepley 
1649371c9d4SSatish Balay static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order) {
1656823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1661e50132fSMatthew G. Knepley 
1671e50132fSMatthew G. Knepley   PetscFunctionBegin;
1686823f3c5SBlaise Bourdin   *order = exo->order;
1696823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
1706823f3c5SBlaise Bourdin }
1716823f3c5SBlaise Bourdin 
1729371c9d4SSatish Balay static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order) {
1736823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1746823f3c5SBlaise Bourdin 
1756823f3c5SBlaise Bourdin   PetscFunctionBegin;
1766823f3c5SBlaise Bourdin   exo->order = order;
1771e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1781e50132fSMatthew G. Knepley }
1791e50132fSMatthew G. Knepley 
1801e50132fSMatthew G. Knepley /*MC
18100373969SVaclav Hapla    PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
1821e50132fSMatthew G. Knepley 
183db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`,
184db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
1851e50132fSMatthew G. Knepley 
1861e50132fSMatthew G. Knepley   Level: beginner
1871e50132fSMatthew G. Knepley M*/
1881e50132fSMatthew G. Knepley 
1899371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v) {
1901e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo;
1911e50132fSMatthew G. Knepley 
1921e50132fSMatthew G. Knepley   PetscFunctionBegin;
193*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&exo));
1941e50132fSMatthew G. Knepley 
1951e50132fSMatthew G. Knepley   v->data                = (void *)exo;
1961e50132fSMatthew G. Knepley   v->ops->destroy        = PetscViewerDestroy_ExodusII;
1971e50132fSMatthew G. Knepley   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
1981e50132fSMatthew G. Knepley   v->ops->setup          = PetscViewerSetUp_ExodusII;
1991e50132fSMatthew G. Knepley   v->ops->view           = PetscViewerView_ExodusII;
2001e50132fSMatthew G. Knepley   v->ops->flush          = 0;
2011e50132fSMatthew G. Knepley   exo->btype             = (PetscFileMode)-1;
2021e50132fSMatthew G. Knepley   exo->filename          = 0;
2031e50132fSMatthew G. Knepley   exo->exoid             = -1;
2041e50132fSMatthew G. Knepley 
2059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_ExodusII));
2069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_ExodusII));
2079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_ExodusII));
2089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_ExodusII));
2099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetId_C", PetscViewerExodusIIGetId_ExodusII));
2109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerSetOrder_C", PetscViewerExodusIISetOrder_ExodusII));
2119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetOrder_C", PetscViewerExodusIIGetOrder_ExodusII));
2121e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
2131e50132fSMatthew G. Knepley }
2141e50132fSMatthew G. Knepley 
2151e50132fSMatthew G. Knepley /*
21678569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
21778569bb4SMatthew G. Knepley 
2186823f3c5SBlaise Bourdin   Collective
21978569bb4SMatthew G. Knepley 
22078569bb4SMatthew G. Knepley   Input Parameters:
22178569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
22278569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
22378569bb4SMatthew G. Knepley - name     - the name of the result
22478569bb4SMatthew G. Knepley 
22578569bb4SMatthew G. Knepley   Output Parameters:
22678569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
22778569bb4SMatthew G. Knepley 
22878569bb4SMatthew G. Knepley   Notes:
22978569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
23078569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
23178569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
23278569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
23378569bb4SMatthew G. Knepley 
23478569bb4SMatthew G. Knepley   Level: beginner
23578569bb4SMatthew G. Knepley 
236c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
23778569bb4SMatthew G. Knepley */
2389371c9d4SSatish Balay PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex) {
2396de834b4SMatthew G. Knepley   int       num_vars, i, j;
24078569bb4SMatthew G. Knepley   char      ext_name[MAX_STR_LENGTH + 1], var_name[MAX_STR_LENGTH + 1];
24178569bb4SMatthew G. Knepley   const int num_suffix = 5;
24278569bb4SMatthew G. Knepley   char     *suffix[5];
24378569bb4SMatthew G. Knepley   PetscBool flg;
24478569bb4SMatthew G. Knepley 
24578569bb4SMatthew G. Knepley   PetscFunctionBegin;
24678569bb4SMatthew G. Knepley   suffix[0] = (char *)"";
24778569bb4SMatthew G. Knepley   suffix[1] = (char *)"_X";
24878569bb4SMatthew G. Knepley   suffix[2] = (char *)"_XX";
24978569bb4SMatthew G. Knepley   suffix[3] = (char *)"_1";
25078569bb4SMatthew G. Knepley   suffix[4] = (char *)"_11";
25178569bb4SMatthew G. Knepley 
2526823f3c5SBlaise Bourdin   *varIndex = -1;
253792fecdfSBarry Smith   PetscCallExternal(ex_get_variable_param, exoid, obj_type, &num_vars);
25478569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
255792fecdfSBarry Smith     PetscCallExternal(ex_get_variable_name, exoid, obj_type, i + 1, var_name);
25678569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j) {
2579566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
2589566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
2599566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(ext_name, var_name, &flg));
260ad540459SPierre Jolivet       if (flg) *varIndex = i + 1;
2616823f3c5SBlaise Bourdin     }
2626823f3c5SBlaise Bourdin   }
26378569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
26478569bb4SMatthew G. Knepley }
26578569bb4SMatthew G. Knepley 
26678569bb4SMatthew G. Knepley /*
2676823f3c5SBlaise Bourdin   DMView_PlexExodusII - Write a DM to disk in exodus format
26878569bb4SMatthew G. Knepley 
26978569bb4SMatthew G. Knepley   Collective on dm
27078569bb4SMatthew G. Knepley 
27178569bb4SMatthew G. Knepley   Input Parameters:
27278569bb4SMatthew G. Knepley + dm  - The dm to be written
2736823f3c5SBlaise Bourdin . viewer - an exodusII viewer
27478569bb4SMatthew G. Knepley 
27578569bb4SMatthew G. Knepley   Notes:
27678569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
27778569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
27878569bb4SMatthew G. Knepley 
2796823f3c5SBlaise Bourdin   If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices)
2806823f3c5SBlaise Bourdin   will be written.
28178569bb4SMatthew G. Knepley 
28278569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
28378569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
2846823f3c5SBlaise Bourdin   The order of the mesh shall be set using PetscViewerExodusIISetOrder
28578569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
28678569bb4SMatthew G. Knepley 
2876823f3c5SBlaise Bourdin   This function will only handle TRI, TET, QUAD, and HEX cells.
28878569bb4SMatthew G. Knepley   Level: beginner
28978569bb4SMatthew G. Knepley 
2906823f3c5SBlaise Bourdin .seealso:
29178569bb4SMatthew G. Knepley */
2929371c9d4SSatish Balay PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer) {
2939371c9d4SSatish Balay   enum ElemType {
2949371c9d4SSatish Balay     TRI,
2959371c9d4SSatish Balay     QUAD,
2969371c9d4SSatish Balay     TET,
2979371c9d4SSatish Balay     HEX
2989371c9d4SSatish Balay   };
29978569bb4SMatthew G. Knepley   MPI_Comm        comm;
3006823f3c5SBlaise Bourdin   PetscInt        degree; /* the order of the mesh */
30178569bb4SMatthew G. Knepley   /* Connectivity Variables */
30278569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
30378569bb4SMatthew G. Knepley   /* Cell Sets */
30478569bb4SMatthew G. Knepley   DMLabel         csLabel;
30578569bb4SMatthew G. Knepley   IS              csIS;
30678569bb4SMatthew G. Knepley   const PetscInt *csIdx;
30778569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
30878569bb4SMatthew G. Knepley   enum ElemType  *type;
309fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
31078569bb4SMatthew G. Knepley   /* Coordinate Variables */
31178569bb4SMatthew G. Knepley   DM              cdm;
3126823f3c5SBlaise Bourdin   PetscSection    coordSection;
31378569bb4SMatthew G. Knepley   Vec             coord;
31478569bb4SMatthew G. Knepley   PetscInt      **nodes;
315eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
31678569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
31778569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
31878569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
31978569bb4SMatthew G. Knepley   const char     *dmName;
320fe209ef9SBlaise Bourdin   PetscInt        nodesTriP1[4]  = {3, 0, 0, 0};
321fe209ef9SBlaise Bourdin   PetscInt        nodesTriP2[4]  = {3, 3, 0, 0};
322fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP1[4] = {4, 0, 0, 0};
323fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP2[4] = {4, 4, 0, 1};
324fe209ef9SBlaise Bourdin   PetscInt        nodesTetP1[4]  = {4, 0, 0, 0};
325fe209ef9SBlaise Bourdin   PetscInt        nodesTetP2[4]  = {4, 6, 0, 0};
326fe209ef9SBlaise Bourdin   PetscInt        nodesHexP1[4]  = {8, 0, 0, 0};
327fe209ef9SBlaise Bourdin   PetscInt        nodesHexP2[4]  = {8, 12, 6, 1};
3286823f3c5SBlaise Bourdin   int             CPU_word_size, IO_word_size, EXO_mode;
3296823f3c5SBlaise Bourdin   float           EXO_version;
3306823f3c5SBlaise Bourdin 
3316823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
3326823f3c5SBlaise Bourdin 
33378569bb4SMatthew G. Knepley   PetscFunctionBegin;
3349566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3376823f3c5SBlaise Bourdin 
3386823f3c5SBlaise Bourdin   /*
3396823f3c5SBlaise Bourdin     Creating coordSection is a collective operation so we do it somewhat out of sequence
3406823f3c5SBlaise Bourdin   */
3419566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &coordSection));
3429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
343c5853193SPierre Jolivet   if (rank == 0) {
3446823f3c5SBlaise Bourdin     switch (exo->btype) {
3456823f3c5SBlaise Bourdin     case FILE_MODE_READ:
3466823f3c5SBlaise Bourdin     case FILE_MODE_APPEND:
3476823f3c5SBlaise Bourdin     case FILE_MODE_UPDATE:
3486823f3c5SBlaise Bourdin     case FILE_MODE_APPEND_UPDATE:
3496823f3c5SBlaise Bourdin       /* exodusII does not allow writing geometry to an existing file */
35098921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
3516823f3c5SBlaise Bourdin     case FILE_MODE_WRITE:
3526823f3c5SBlaise Bourdin       /* Create an empty file if one already exists*/
3536823f3c5SBlaise Bourdin       EXO_mode = EX_CLOBBER;
3546823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
3556823f3c5SBlaise Bourdin       EXO_mode += EX_ALL_INT64_API;
3566823f3c5SBlaise Bourdin #endif
3576823f3c5SBlaise Bourdin       CPU_word_size = sizeof(PetscReal);
3586823f3c5SBlaise Bourdin       IO_word_size  = sizeof(PetscReal);
3596823f3c5SBlaise Bourdin       exo->exoid    = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
36008401ef6SPierre Jolivet       PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
3616823f3c5SBlaise Bourdin 
3626823f3c5SBlaise Bourdin       break;
3639371c9d4SSatish Balay     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3646823f3c5SBlaise Bourdin     }
3656823f3c5SBlaise Bourdin 
36678569bb4SMatthew G. Knepley     /* --- Get DM info --- */
3679566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)dm, &dmName));
3689566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepth(dm, &depth));
3699566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
3709566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3719566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
3739566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
3749566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
37578569bb4SMatthew G. Knepley     numCells    = cEnd - cStart;
37678569bb4SMatthew G. Knepley     numEdges    = eEnd - eStart;
37778569bb4SMatthew G. Knepley     numVertices = vEnd - vStart;
3789371c9d4SSatish Balay     if (depth == 3) {
3799371c9d4SSatish Balay       numFaces = fEnd - fStart;
3809371c9d4SSatish Balay     } else {
3819371c9d4SSatish Balay       numFaces = 0;
3829371c9d4SSatish Balay     }
3839566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs));
3849566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs));
3859566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs));
3869566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coord));
3879566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
38878569bb4SMatthew G. Knepley     if (num_cs > 0) {
3899566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel));
3909566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(csLabel, &csIS));
3919566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(csIS, &csIdx));
39278569bb4SMatthew G. Knepley     }
3939566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &nodes));
39478569bb4SMatthew G. Knepley     /* Set element type for each block and compute total number of nodes */
3959566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &type));
39678569bb4SMatthew G. Knepley     numNodes = numVertices;
3976823f3c5SBlaise Bourdin 
3989566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIGetOrder(viewer, &degree));
399ad540459SPierre Jolivet     if (degree == 2) numNodes += numEdges;
40078569bb4SMatthew G. Knepley     cellsNotInConnectivity = numCells;
40178569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
40278569bb4SMatthew G. Knepley       IS              stratumIS;
40378569bb4SMatthew G. Knepley       const PetscInt *cells;
40478569bb4SMatthew G. Knepley       PetscScalar    *xyz = NULL;
40578569bb4SMatthew G. Knepley       PetscInt        csSize, closureSize;
40678569bb4SMatthew G. Knepley 
4079566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4089566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4099566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
4109566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
41178569bb4SMatthew G. Knepley       switch (dim) {
41278569bb4SMatthew G. Knepley       case 2:
4139371c9d4SSatish Balay         if (closureSize == 3 * dim) {
4149371c9d4SSatish Balay           type[cs] = TRI;
4159371c9d4SSatish Balay         } else if (closureSize == 4 * dim) {
4169371c9d4SSatish Balay           type[cs] = QUAD;
4179371c9d4SSatish Balay         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim);
41878569bb4SMatthew G. Knepley         break;
41978569bb4SMatthew G. Knepley       case 3:
4209371c9d4SSatish Balay         if (closureSize == 4 * dim) {
4219371c9d4SSatish Balay           type[cs] = TET;
4229371c9d4SSatish Balay         } else if (closureSize == 8 * dim) {
4239371c9d4SSatish Balay           type[cs] = HEX;
4249371c9d4SSatish Balay         } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %" PetscInt_FMT " in dimension %" PetscInt_FMT " has no ExodusII type", closureSize / dim, dim);
42578569bb4SMatthew G. Knepley         break;
42663a3b9bcSJacob Faibussowitsch       default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim);
42778569bb4SMatthew G. Knepley       }
428ad540459SPierre Jolivet       if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize;
4299371c9d4SSatish Balay       if ((degree == 2) && (type[cs] == HEX)) {
4309371c9d4SSatish Balay         numNodes += csSize;
4319371c9d4SSatish Balay         numNodes += numFaces;
4329371c9d4SSatish Balay       }
4339566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
43478569bb4SMatthew G. Knepley       /* Set nodes and Element type */
43578569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
43678569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesTriP1;
43778569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTriP2;
43878569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
43978569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesQuadP1;
44078569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesQuadP2;
44178569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
44278569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesTetP1;
44378569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTetP2;
44478569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
44578569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesHexP1;
44678569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesHexP2;
44778569bb4SMatthew G. Knepley       }
44878569bb4SMatthew G. Knepley       /* Compute the number of cells not in the connectivity table */
44978569bb4SMatthew G. Knepley       cellsNotInConnectivity -= nodes[cs][3] * csSize;
45078569bb4SMatthew G. Knepley 
4519566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
4529566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
45378569bb4SMatthew G. Knepley     }
45448a46eb9SPierre Jolivet     if (num_cs > 0) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);
45578569bb4SMatthew G. Knepley     /* --- Connectivity --- */
45678569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
45778569bb4SMatthew G. Knepley       IS              stratumIS;
45878569bb4SMatthew G. Knepley       const PetscInt *cells;
459fe209ef9SBlaise Bourdin       PetscInt       *connect, off = 0;
460eae8a387SMatthew G. Knepley       PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
46178569bb4SMatthew G. Knepley       PetscInt        csSize, c, connectSize, closureSize;
46278569bb4SMatthew G. Knepley       char           *elem_type        = NULL;
46378569bb4SMatthew G. Knepley       char            elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4";
46478569bb4SMatthew G. Knepley       char            elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9";
46578569bb4SMatthew G. Knepley       char            elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8";
46678569bb4SMatthew G. Knepley       char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
46778569bb4SMatthew G. Knepley 
4689566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4699566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4709566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
47178569bb4SMatthew G. Knepley       /* Set Element type */
47278569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
47378569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_tri3;
47478569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tri6;
47578569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
47678569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_quad4;
47778569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_quad9;
47878569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
47978569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_tet4;
48078569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tet10;
48178569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
48278569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_hex8;
48378569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_hex27;
48478569bb4SMatthew G. Knepley       }
48578569bb4SMatthew G. Knepley       connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
4869566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect));
487792fecdfSBarry Smith       PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1);
48878569bb4SMatthew G. Knepley       /* Find number of vertices, edges, and faces in the closure */
48978569bb4SMatthew G. Knepley       verticesInClosure = nodes[cs][0];
49078569bb4SMatthew G. Knepley       if (depth > 1) {
49178569bb4SMatthew G. Knepley         if (dim == 2) {
4929566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure));
49378569bb4SMatthew G. Knepley         } else if (dim == 3) {
49478569bb4SMatthew G. Knepley           PetscInt *closure = NULL;
49578569bb4SMatthew G. Knepley 
4969566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure));
4979566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
49878569bb4SMatthew G. Knepley           edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
4999566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
50078569bb4SMatthew G. Knepley         }
50178569bb4SMatthew G. Knepley       }
50278569bb4SMatthew G. Knepley       /* Get connectivity for each cell */
50378569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
50478569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
50578569bb4SMatthew G. Knepley         PetscInt  temp, i;
50678569bb4SMatthew G. Knepley 
5079566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
50878569bb4SMatthew G. Knepley         for (i = 0; i < connectSize; ++i) {
50978569bb4SMatthew G. Knepley           if (i < nodes[cs][0]) { /* Vertices */
510fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1;
511fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
51278569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */
513fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1;
514fe209ef9SBlaise Bourdin             if (nodes[cs][2] == 0) connect[i + off] -= numFaces;
515fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
51678569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */
517fe209ef9SBlaise Bourdin             connect[i + off] = closure[0] + 1;
518fe209ef9SBlaise Bourdin             connect[i + off] -= skipCells;
51978569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */
520fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1;
521fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
52278569bb4SMatthew G. Knepley           } else {
523fe209ef9SBlaise Bourdin             connect[i + off] = -1;
52478569bb4SMatthew G. Knepley           }
52578569bb4SMatthew G. Knepley         }
52678569bb4SMatthew G. Knepley         /* Tetrahedra are inverted */
52778569bb4SMatthew G. Knepley         if (type[cs] == TET) {
5289371c9d4SSatish Balay           temp             = connect[0 + off];
5299371c9d4SSatish Balay           connect[0 + off] = connect[1 + off];
5309371c9d4SSatish Balay           connect[1 + off] = temp;
53178569bb4SMatthew G. Knepley           if (degree == 2) {
5329371c9d4SSatish Balay             temp             = connect[5 + off];
5339371c9d4SSatish Balay             connect[5 + off] = connect[6 + off];
5349371c9d4SSatish Balay             connect[6 + off] = temp;
5359371c9d4SSatish Balay             temp             = connect[7 + off];
5369371c9d4SSatish Balay             connect[7 + off] = connect[8 + off];
5379371c9d4SSatish Balay             connect[8 + off] = temp;
53878569bb4SMatthew G. Knepley           }
53978569bb4SMatthew G. Knepley         }
54078569bb4SMatthew G. Knepley         /* Hexahedra are inverted */
54178569bb4SMatthew G. Knepley         if (type[cs] == HEX) {
5429371c9d4SSatish Balay           temp             = connect[1 + off];
5439371c9d4SSatish Balay           connect[1 + off] = connect[3 + off];
5449371c9d4SSatish Balay           connect[3 + off] = temp;
54578569bb4SMatthew G. Knepley           if (degree == 2) {
5469371c9d4SSatish Balay             temp              = connect[8 + off];
5479371c9d4SSatish Balay             connect[8 + off]  = connect[11 + off];
5489371c9d4SSatish Balay             connect[11 + off] = temp;
5499371c9d4SSatish Balay             temp              = connect[9 + off];
5509371c9d4SSatish Balay             connect[9 + off]  = connect[10 + off];
5519371c9d4SSatish Balay             connect[10 + off] = temp;
5529371c9d4SSatish Balay             temp              = connect[16 + off];
5539371c9d4SSatish Balay             connect[16 + off] = connect[17 + off];
5549371c9d4SSatish Balay             connect[17 + off] = temp;
5559371c9d4SSatish Balay             temp              = connect[18 + off];
5569371c9d4SSatish Balay             connect[18 + off] = connect[19 + off];
5579371c9d4SSatish Balay             connect[19 + off] = temp;
55878569bb4SMatthew G. Knepley 
5599371c9d4SSatish Balay             temp              = connect[12 + off];
5609371c9d4SSatish Balay             connect[12 + off] = connect[16 + off];
5619371c9d4SSatish Balay             connect[16 + off] = temp;
5629371c9d4SSatish Balay             temp              = connect[13 + off];
5639371c9d4SSatish Balay             connect[13 + off] = connect[17 + off];
5649371c9d4SSatish Balay             connect[17 + off] = temp;
5659371c9d4SSatish Balay             temp              = connect[14 + off];
5669371c9d4SSatish Balay             connect[14 + off] = connect[18 + off];
5679371c9d4SSatish Balay             connect[18 + off] = temp;
5689371c9d4SSatish Balay             temp              = connect[15 + off];
5699371c9d4SSatish Balay             connect[15 + off] = connect[19 + off];
5709371c9d4SSatish Balay             connect[19 + off] = temp;
57178569bb4SMatthew G. Knepley 
5729371c9d4SSatish Balay             temp              = connect[23 + off];
5739371c9d4SSatish Balay             connect[23 + off] = connect[26 + off];
5749371c9d4SSatish Balay             connect[26 + off] = temp;
5759371c9d4SSatish Balay             temp              = connect[24 + off];
5769371c9d4SSatish Balay             connect[24 + off] = connect[25 + off];
5779371c9d4SSatish Balay             connect[25 + off] = temp;
5789371c9d4SSatish Balay             temp              = connect[25 + off];
5799371c9d4SSatish Balay             connect[25 + off] = connect[26 + off];
5809371c9d4SSatish Balay             connect[26 + off] = temp;
58178569bb4SMatthew G. Knepley           }
58278569bb4SMatthew G. Knepley         }
583fe209ef9SBlaise Bourdin         off += connectSize;
5849566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
58578569bb4SMatthew G. Knepley       }
586792fecdfSBarry Smith       PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0);
58778569bb4SMatthew G. Knepley       skipCells += (nodes[cs][3] == 0) * csSize;
5889566063dSJacob Faibussowitsch       PetscCall(PetscFree(connect));
5899566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
5909566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
59178569bb4SMatthew G. Knepley     }
5929566063dSJacob Faibussowitsch     PetscCall(PetscFree(type));
59378569bb4SMatthew G. Knepley     /* --- Coordinates --- */
5949566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd));
5951e50132fSMatthew G. Knepley     if (num_cs) {
59678569bb4SMatthew G. Knepley       for (d = 0; d < depth; ++d) {
5979566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
59848a46eb9SPierre Jolivet         for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0));
59978569bb4SMatthew G. Knepley       }
6001e50132fSMatthew G. Knepley     }
60178569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
60278569bb4SMatthew G. Knepley       IS              stratumIS;
60378569bb4SMatthew G. Knepley       const PetscInt *cells;
60478569bb4SMatthew G. Knepley       PetscInt        csSize, c;
60578569bb4SMatthew G. Knepley 
6069566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
6079566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
6089566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
60948a46eb9SPierre Jolivet       for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0));
6109566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
6119566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
61278569bb4SMatthew G. Knepley     }
61378569bb4SMatthew G. Knepley     if (num_cs > 0) {
6149566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(csIS, &csIdx));
6159566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&csIS));
61678569bb4SMatthew G. Knepley     }
6179566063dSJacob Faibussowitsch     PetscCall(PetscFree(nodes));
6189566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(coordSection));
61978569bb4SMatthew G. Knepley     if (numNodes > 0) {
62078569bb4SMatthew G. Knepley       const char  *coordNames[3] = {"x", "y", "z"};
621233c95e0SFrancesco Ballarin       PetscScalar *closure, *cval;
622233c95e0SFrancesco Ballarin       PetscReal   *coords;
62378569bb4SMatthew G. Knepley       PetscInt     hasDof, n = 0;
62478569bb4SMatthew G. Knepley 
6256823f3c5SBlaise Bourdin       /* There can't be more than 24 values in the closure of a point for the coord coordSection */
6269566063dSJacob Faibussowitsch       PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure));
6279566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord));
6289566063dSJacob Faibussowitsch       PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
62978569bb4SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
6309566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(coordSection, p, &hasDof));
63178569bb4SMatthew G. Knepley         if (hasDof) {
63278569bb4SMatthew G. Knepley           PetscInt closureSize = 24, j;
63378569bb4SMatthew G. Knepley 
6349566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure));
63578569bb4SMatthew G. Knepley           for (d = 0; d < dim; ++d) {
63678569bb4SMatthew G. Knepley             cval[d] = 0.0;
63778569bb4SMatthew G. Knepley             for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d];
638233c95e0SFrancesco Ballarin             coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize;
63978569bb4SMatthew G. Knepley           }
64078569bb4SMatthew G. Knepley           ++n;
64178569bb4SMatthew G. Knepley         }
64278569bb4SMatthew G. Knepley       }
643792fecdfSBarry Smith       PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]);
6449566063dSJacob Faibussowitsch       PetscCall(PetscFree3(coords, cval, closure));
645792fecdfSBarry Smith       PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames);
64678569bb4SMatthew G. Knepley     }
6476823f3c5SBlaise Bourdin 
648fe209ef9SBlaise Bourdin     /* --- Node Sets/Vertex Sets --- */
649fe209ef9SBlaise Bourdin     DMHasLabel(dm, "Vertex Sets", &hasLabel);
650fe209ef9SBlaise Bourdin     if (hasLabel) {
651fe209ef9SBlaise Bourdin       PetscInt        i, vs, vsSize;
652fe209ef9SBlaise Bourdin       const PetscInt *vsIdx, *vertices;
653fe209ef9SBlaise Bourdin       PetscInt       *nodeList;
654fe209ef9SBlaise Bourdin       IS              vsIS, stratumIS;
655fe209ef9SBlaise Bourdin       DMLabel         vsLabel;
6569566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel));
6579566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(vsLabel, &vsIS));
6589566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(vsIS, &vsIdx));
659fe209ef9SBlaise Bourdin       for (vs = 0; vs < num_vs; ++vs) {
6609566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS));
6619566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(stratumIS, &vertices));
6629566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &vsSize));
6639566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(vsSize, &nodeList));
664ad540459SPierre Jolivet         for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1;
665792fecdfSBarry Smith         PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0);
666792fecdfSBarry Smith         PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL);
6679566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(stratumIS, &vertices));
6689566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
6699566063dSJacob Faibussowitsch         PetscCall(PetscFree(nodeList));
670fe209ef9SBlaise Bourdin       }
6719566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(vsIS, &vsIdx));
6729566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&vsIS));
673fe209ef9SBlaise Bourdin     }
674fe209ef9SBlaise Bourdin     /* --- Side Sets/Face Sets --- */
6759566063dSJacob Faibussowitsch     PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel));
676fe209ef9SBlaise Bourdin     if (hasLabel) {
677fe209ef9SBlaise Bourdin       PetscInt        i, j, fs, fsSize;
678fe209ef9SBlaise Bourdin       const PetscInt *fsIdx, *faces;
679fe209ef9SBlaise Bourdin       IS              fsIS, stratumIS;
680fe209ef9SBlaise Bourdin       DMLabel         fsLabel;
681fe209ef9SBlaise Bourdin       PetscInt        numPoints, *points;
682fe209ef9SBlaise Bourdin       PetscInt        elem_list_size = 0;
683fe209ef9SBlaise Bourdin       PetscInt       *elem_list, *elem_ind, *side_list;
684fe209ef9SBlaise Bourdin 
6859566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel));
686fe209ef9SBlaise Bourdin       /* Compute size of Node List and Element List */
6879566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(fsLabel, &fsIS));
6889566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fsIS, &fsIdx));
689fe209ef9SBlaise Bourdin       for (fs = 0; fs < num_fs; ++fs) {
6909566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
6919566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &fsSize));
692fe209ef9SBlaise Bourdin         elem_list_size += fsSize;
6939566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
694fe209ef9SBlaise Bourdin       }
695d0609cedSBarry Smith       PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list));
696fe209ef9SBlaise Bourdin       elem_ind[0] = 0;
697fe209ef9SBlaise Bourdin       for (fs = 0; fs < num_fs; ++fs) {
6989566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
6999566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(stratumIS, &faces));
7009566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &fsSize));
701fe209ef9SBlaise Bourdin         /* Set Parameters */
702792fecdfSBarry Smith         PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0);
703fe209ef9SBlaise Bourdin         /* Indices */
704ad540459SPierre Jolivet         if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize;
705fe209ef9SBlaise Bourdin 
706fe209ef9SBlaise Bourdin         for (i = 0; i < fsSize; ++i) {
707fe209ef9SBlaise Bourdin           /* Element List */
708fe209ef9SBlaise Bourdin           points = NULL;
7099566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
710fe209ef9SBlaise Bourdin           elem_list[elem_ind[fs] + i] = points[2] + 1;
7119566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
712fe209ef9SBlaise Bourdin 
713fe209ef9SBlaise Bourdin           /* Side List */
714fe209ef9SBlaise Bourdin           points = NULL;
7159566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
716fe209ef9SBlaise Bourdin           for (j = 1; j < numPoints; ++j) {
717ad540459SPierre Jolivet             if (points[j * 2] == faces[i]) break;
718fe209ef9SBlaise Bourdin           }
719fe209ef9SBlaise Bourdin           /* Convert HEX sides */
720fe209ef9SBlaise Bourdin           if (numPoints == 27) {
7219371c9d4SSatish Balay             if (j == 1) {
7229371c9d4SSatish Balay               j = 5;
7239371c9d4SSatish Balay             } else if (j == 2) {
7249371c9d4SSatish Balay               j = 6;
7259371c9d4SSatish Balay             } else if (j == 3) {
7269371c9d4SSatish Balay               j = 1;
7279371c9d4SSatish Balay             } else if (j == 4) {
7289371c9d4SSatish Balay               j = 3;
7299371c9d4SSatish Balay             } else if (j == 5) {
7309371c9d4SSatish Balay               j = 2;
7319371c9d4SSatish Balay             } else if (j == 6) {
7329371c9d4SSatish Balay               j = 4;
7339371c9d4SSatish Balay             }
734fe209ef9SBlaise Bourdin           }
735fe209ef9SBlaise Bourdin           /* Convert TET sides */
736fe209ef9SBlaise Bourdin           if (numPoints == 15) {
737fe209ef9SBlaise Bourdin             --j;
738ad540459SPierre Jolivet             if (j == 0) j = 4;
739fe209ef9SBlaise Bourdin           }
740fe209ef9SBlaise Bourdin           side_list[elem_ind[fs] + i] = j;
7419566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
742fe209ef9SBlaise Bourdin         }
7439566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(stratumIS, &faces));
7449566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
745fe209ef9SBlaise Bourdin       }
7469566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fsIS, &fsIdx));
7479566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&fsIS));
748fe209ef9SBlaise Bourdin 
749fe209ef9SBlaise Bourdin       /* Put side sets */
75048a46eb9SPierre Jolivet       for (fs = 0; fs < num_fs; ++fs) PetscCallExternal(ex_put_set, exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]);
7519566063dSJacob Faibussowitsch       PetscCall(PetscFree3(elem_ind, elem_list, side_list));
752fe209ef9SBlaise Bourdin     }
7536823f3c5SBlaise Bourdin     /*
7546823f3c5SBlaise Bourdin       close the exodus file
7556823f3c5SBlaise Bourdin     */
7566823f3c5SBlaise Bourdin     ex_close(exo->exoid);
7576823f3c5SBlaise Bourdin     exo->exoid = -1;
7586823f3c5SBlaise Bourdin   }
7599566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&coordSection));
7606823f3c5SBlaise Bourdin 
7616823f3c5SBlaise Bourdin   /*
7626823f3c5SBlaise Bourdin     reopen the file in parallel
7636823f3c5SBlaise Bourdin   */
7646823f3c5SBlaise Bourdin   EXO_mode = EX_WRITE;
7656823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
7666823f3c5SBlaise Bourdin   EXO_mode += EX_ALL_INT64_API;
7676823f3c5SBlaise Bourdin #endif
7686823f3c5SBlaise Bourdin   CPU_word_size = sizeof(PetscReal);
7696823f3c5SBlaise Bourdin   IO_word_size  = sizeof(PetscReal);
7706823f3c5SBlaise Bourdin   exo->exoid    = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL);
77108401ef6SPierre Jolivet   PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
7726823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
7736823f3c5SBlaise Bourdin }
7746823f3c5SBlaise Bourdin 
7756823f3c5SBlaise Bourdin /*
7766823f3c5SBlaise Bourdin   VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
7776823f3c5SBlaise Bourdin 
7786823f3c5SBlaise Bourdin   Collective on v
7796823f3c5SBlaise Bourdin 
7806823f3c5SBlaise Bourdin   Input Parameters:
7816823f3c5SBlaise Bourdin + v  - The vector to be written
7826823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
7836823f3c5SBlaise Bourdin 
7846823f3c5SBlaise Bourdin   Notes:
7856823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
7866823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
7876823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
7886823f3c5SBlaise Bourdin   amongst all variables.
7896823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
7906823f3c5SBlaise Bourdin   possibly corrupting the file
7916823f3c5SBlaise Bourdin 
7926823f3c5SBlaise Bourdin   Level: beginner
7936823f3c5SBlaise Bourdin 
794c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
7956823f3c5SBlaise Bourdin @*/
7969371c9d4SSatish Balay PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer) {
7976823f3c5SBlaise Bourdin   DM          dm;
7986823f3c5SBlaise Bourdin   MPI_Comm    comm;
7996823f3c5SBlaise Bourdin   PetscMPIInt rank;
8006823f3c5SBlaise Bourdin 
8016823f3c5SBlaise Bourdin   int         exoid, offsetN = 0, offsetZ = 0;
8026823f3c5SBlaise Bourdin   const char *vecname;
8036823f3c5SBlaise Bourdin   PetscInt    step;
8046823f3c5SBlaise Bourdin 
8056823f3c5SBlaise Bourdin   PetscFunctionBegin;
8069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
8079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8089566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
8099566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8109566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
8116823f3c5SBlaise Bourdin 
8129566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
8139566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN));
8149566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ));
8151dca8a05SBarry Smith   PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
8166823f3c5SBlaise Bourdin   if (offsetN > 0) {
8179566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN));
8186823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
8199566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ));
82098921bdaSJacob Faibussowitsch   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
8216823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
8226823f3c5SBlaise Bourdin }
8236823f3c5SBlaise Bourdin 
8246823f3c5SBlaise Bourdin /*
8256823f3c5SBlaise Bourdin   VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
8266823f3c5SBlaise Bourdin 
8276823f3c5SBlaise Bourdin   Collective on v
8286823f3c5SBlaise Bourdin 
8296823f3c5SBlaise Bourdin   Input Parameters:
8306823f3c5SBlaise Bourdin + v  - The vector to be written
8316823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
8326823f3c5SBlaise Bourdin 
8336823f3c5SBlaise Bourdin   Notes:
8346823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
8356823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
8366823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
8376823f3c5SBlaise Bourdin   amongst all variables.
8386823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
8396823f3c5SBlaise Bourdin   possibly corrupting the file
8406823f3c5SBlaise Bourdin 
8416823f3c5SBlaise Bourdin   Level: beginner
8426823f3c5SBlaise Bourdin 
843c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
8446823f3c5SBlaise Bourdin @*/
8459371c9d4SSatish Balay PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer) {
8466823f3c5SBlaise Bourdin   DM          dm;
8476823f3c5SBlaise Bourdin   MPI_Comm    comm;
8486823f3c5SBlaise Bourdin   PetscMPIInt rank;
8496823f3c5SBlaise Bourdin 
8506823f3c5SBlaise Bourdin   int         exoid, offsetN = 0, offsetZ = 0;
8516823f3c5SBlaise Bourdin   const char *vecname;
8526823f3c5SBlaise Bourdin   PetscInt    step;
8536823f3c5SBlaise Bourdin 
8546823f3c5SBlaise Bourdin   PetscFunctionBegin;
8559566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
8569566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8579566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
8589566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
8606823f3c5SBlaise Bourdin 
8619566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
8629566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN));
8639566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ));
8641dca8a05SBarry Smith   PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
8651dca8a05SBarry Smith   if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN));
8661dca8a05SBarry Smith   else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ));
8671dca8a05SBarry Smith   else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
86878569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
86978569bb4SMatthew G. Knepley }
87078569bb4SMatthew G. Knepley 
87178569bb4SMatthew G. Knepley /*
87278569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
87378569bb4SMatthew G. Knepley 
87478569bb4SMatthew G. Knepley   Collective on v
87578569bb4SMatthew G. Knepley 
87678569bb4SMatthew G. Knepley   Input Parameters:
87778569bb4SMatthew G. Knepley + v  - The vector to be written
87878569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
8796823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
8806823f3c5SBlaise Bourdin - offset - the location of the variable in the file
88178569bb4SMatthew G. Knepley 
88278569bb4SMatthew G. Knepley   Notes:
88378569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
88478569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
88578569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
88678569bb4SMatthew G. Knepley   amongst all nodal variables.
88778569bb4SMatthew G. Knepley 
88878569bb4SMatthew G. Knepley   Level: beginner
88978569bb4SMatthew G. Knepley 
890c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
89178569bb4SMatthew G. Knepley @*/
8929371c9d4SSatish Balay PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) {
89378569bb4SMatthew G. Knepley   MPI_Comm           comm;
89478569bb4SMatthew G. Knepley   PetscMPIInt        size;
89578569bb4SMatthew G. Knepley   DM                 dm;
89678569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
89722a7544eSVaclav Hapla   const PetscScalar *varray;
89878569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
89978569bb4SMatthew G. Knepley   PetscBool          useNatural;
90078569bb4SMatthew G. Knepley 
90178569bb4SMatthew G. Knepley   PetscFunctionBegin;
9029566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
9039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9049566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
9059566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
90678569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
90778569bb4SMatthew G. Knepley   if (useNatural) {
908b7352c5cSAlexis Marboeuf     PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
9099566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
9109566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
91178569bb4SMatthew G. Knepley   } else {
91278569bb4SMatthew G. Knepley     vNatural = v;
91378569bb4SMatthew G. Knepley   }
9146823f3c5SBlaise Bourdin 
91578569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
91678569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
91778569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
9189566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
9199566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
92078569bb4SMatthew G. Knepley   if (bs == 1) {
9219566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(vNatural, &varray));
922792fecdfSBarry Smith     PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
9239566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(vNatural, &varray));
92478569bb4SMatthew G. Knepley   } else {
92578569bb4SMatthew G. Knepley     IS       compIS;
92678569bb4SMatthew G. Knepley     PetscInt c;
92778569bb4SMatthew G. Knepley 
9289566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
92978569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
9309566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
9319566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
9329566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vComp, &varray));
933792fecdfSBarry Smith       PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
9349566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vComp, &varray));
9359566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
93678569bb4SMatthew G. Knepley     }
9379566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
93878569bb4SMatthew G. Knepley   }
939b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(VecDestroy(&vNatural));
94078569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
94178569bb4SMatthew G. Knepley }
94278569bb4SMatthew G. Knepley 
94378569bb4SMatthew G. Knepley /*
94478569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
94578569bb4SMatthew G. Knepley 
94678569bb4SMatthew G. Knepley   Collective on v
94778569bb4SMatthew G. Knepley 
94878569bb4SMatthew G. Knepley   Input Parameters:
94978569bb4SMatthew G. Knepley + v  - The vector to be written
95078569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
9516823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
9526823f3c5SBlaise Bourdin - offset - the location of the variable in the file
95378569bb4SMatthew G. Knepley 
95478569bb4SMatthew G. Knepley   Notes:
95578569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
95678569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
95778569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
95878569bb4SMatthew G. Knepley   amongst all nodal variables.
95978569bb4SMatthew G. Knepley 
96078569bb4SMatthew G. Knepley   Level: beginner
96178569bb4SMatthew G. Knepley 
962db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
96378569bb4SMatthew G. Knepley */
9649371c9d4SSatish Balay PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset) {
96578569bb4SMatthew G. Knepley   MPI_Comm     comm;
96678569bb4SMatthew G. Knepley   PetscMPIInt  size;
96778569bb4SMatthew G. Knepley   DM           dm;
96878569bb4SMatthew G. Knepley   Vec          vNatural, vComp;
96922a7544eSVaclav Hapla   PetscScalar *varray;
97078569bb4SMatthew G. Knepley   PetscInt     xs, xe, bs;
97178569bb4SMatthew G. Knepley   PetscBool    useNatural;
97278569bb4SMatthew G. Knepley 
97378569bb4SMatthew G. Knepley   PetscFunctionBegin;
9749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
9759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9769566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
9779566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
97878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
979b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
980ad540459SPierre Jolivet   else vNatural = v;
9816823f3c5SBlaise Bourdin 
98278569bb4SMatthew G. Knepley   /* Read local chunk from the file */
9839566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
9849566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
98578569bb4SMatthew G. Knepley   if (bs == 1) {
9869566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vNatural, &varray));
987792fecdfSBarry Smith     PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
9889566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vNatural, &varray));
98978569bb4SMatthew G. Knepley   } else {
99078569bb4SMatthew G. Knepley     IS       compIS;
99178569bb4SMatthew G. Knepley     PetscInt c;
99278569bb4SMatthew G. Knepley 
9939566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
99478569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
9959566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
9969566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
9979566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vComp, &varray));
998792fecdfSBarry Smith       PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
9999566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vComp, &varray));
10009566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
100178569bb4SMatthew G. Knepley     }
10029566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
100378569bb4SMatthew G. Knepley   }
100478569bb4SMatthew G. Knepley   if (useNatural) {
10059566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
10069566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
1007b7352c5cSAlexis Marboeuf     PetscCall(VecDestroy(&vNatural));
100878569bb4SMatthew G. Knepley   }
100978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
101078569bb4SMatthew G. Knepley }
101178569bb4SMatthew G. Knepley 
101278569bb4SMatthew G. Knepley /*
101378569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
101478569bb4SMatthew G. Knepley 
101578569bb4SMatthew G. Knepley   Collective on v
101678569bb4SMatthew G. Knepley 
101778569bb4SMatthew G. Knepley   Input Parameters:
101878569bb4SMatthew G. Knepley + v  - The vector to be written
101978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
10206823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
10216823f3c5SBlaise Bourdin - offset - the location of the variable in the file
102278569bb4SMatthew G. Knepley 
102378569bb4SMatthew G. Knepley   Notes:
102478569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
102578569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
102678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
102778569bb4SMatthew G. Knepley   amongst all zonal variables.
102878569bb4SMatthew G. Knepley 
102978569bb4SMatthew G. Knepley   Level: beginner
103078569bb4SMatthew G. Knepley 
1031c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
103278569bb4SMatthew G. Knepley */
10339371c9d4SSatish Balay PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) {
103478569bb4SMatthew G. Knepley   MPI_Comm           comm;
103578569bb4SMatthew G. Knepley   PetscMPIInt        size;
103678569bb4SMatthew G. Knepley   DM                 dm;
103778569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
103822a7544eSVaclav Hapla   const PetscScalar *varray;
103978569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
104078569bb4SMatthew G. Knepley   PetscBool          useNatural;
104178569bb4SMatthew G. Knepley   IS                 compIS;
104278569bb4SMatthew G. Knepley   PetscInt          *csSize, *csID;
104378569bb4SMatthew G. Knepley   PetscInt           numCS, set, csxs = 0;
104478569bb4SMatthew G. Knepley 
104578569bb4SMatthew G. Knepley   PetscFunctionBegin;
10469566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
10479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10489566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
10499566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
105078569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
105178569bb4SMatthew G. Knepley   if (useNatural) {
1052b7352c5cSAlexis Marboeuf     PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
10539566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
10549566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
105578569bb4SMatthew G. Knepley   } else {
105678569bb4SMatthew G. Knepley     vNatural = v;
105778569bb4SMatthew G. Knepley   }
10586823f3c5SBlaise Bourdin 
105978569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
106078569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
106178569bb4SMatthew G. Knepley      We assume that they are stored sequentially
106278569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1063a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
106478569bb4SMatthew G. Knepley      to figure out what to save where. */
106578569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
10669566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1067792fecdfSBarry Smith   PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
106878569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
106978569bb4SMatthew G. Knepley     ex_block block;
107078569bb4SMatthew G. Knepley 
107178569bb4SMatthew G. Knepley     block.id   = csID[set];
107278569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1073792fecdfSBarry Smith     PetscCallExternal(ex_get_block_param, exoid, &block);
107478569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
107578569bb4SMatthew G. Knepley   }
10769566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
10779566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
10789566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
107978569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
108078569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
108178569bb4SMatthew G. Knepley 
108278569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
108378569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
108478569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
108578569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
108678569bb4SMatthew G. Knepley     if (bs == 1) {
10879566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vNatural, &varray));
1088792fecdfSBarry Smith       PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]);
10899566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vNatural, &varray));
109078569bb4SMatthew G. Knepley     } else {
109178569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
10929566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
10939566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
10949566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vComp, &varray));
1095792fecdfSBarry Smith         PetscCallExternal(ex_put_partial_var, exoid, step, EX_ELEM_BLOCK, offset + c, csID[set], PetscMax(xs / bs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs / bs)]);
10969566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vComp, &varray));
10979566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
109878569bb4SMatthew G. Knepley       }
109978569bb4SMatthew G. Knepley     }
110078569bb4SMatthew G. Knepley     csxs += csSize[set];
110178569bb4SMatthew G. Knepley   }
11029566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
11039566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
1104b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(VecDestroy(&vNatural));
110578569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
110678569bb4SMatthew G. Knepley }
110778569bb4SMatthew G. Knepley 
110878569bb4SMatthew G. Knepley /*
110978569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
111078569bb4SMatthew G. Knepley 
111178569bb4SMatthew G. Knepley   Collective on v
111278569bb4SMatthew G. Knepley 
111378569bb4SMatthew G. Knepley   Input Parameters:
111478569bb4SMatthew G. Knepley + v  - The vector to be written
111578569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
11166823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
11176823f3c5SBlaise Bourdin - offset - the location of the variable in the file
111878569bb4SMatthew G. Knepley 
111978569bb4SMatthew G. Knepley   Notes:
112078569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
112178569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
112278569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
112378569bb4SMatthew G. Knepley   amongst all zonal variables.
112478569bb4SMatthew G. Knepley 
112578569bb4SMatthew G. Knepley   Level: beginner
112678569bb4SMatthew G. Knepley 
1127db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
112878569bb4SMatthew G. Knepley */
11299371c9d4SSatish Balay PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset) {
113078569bb4SMatthew G. Knepley   MPI_Comm     comm;
113178569bb4SMatthew G. Knepley   PetscMPIInt  size;
113278569bb4SMatthew G. Knepley   DM           dm;
113378569bb4SMatthew G. Knepley   Vec          vNatural, vComp;
113422a7544eSVaclav Hapla   PetscScalar *varray;
113578569bb4SMatthew G. Knepley   PetscInt     xs, xe, bs;
113678569bb4SMatthew G. Knepley   PetscBool    useNatural;
113778569bb4SMatthew G. Knepley   IS           compIS;
113878569bb4SMatthew G. Knepley   PetscInt    *csSize, *csID;
113978569bb4SMatthew G. Knepley   PetscInt     numCS, set, csxs = 0;
114078569bb4SMatthew G. Knepley 
114178569bb4SMatthew G. Knepley   PetscFunctionBegin;
11429566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
11439566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
11449566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
11459566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
114678569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
1147b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
1148ad540459SPierre Jolivet   else vNatural = v;
11496823f3c5SBlaise Bourdin 
115078569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
115178569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
115278569bb4SMatthew G. Knepley      We assume that they are stored sequentially
115378569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1154a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
115578569bb4SMatthew G. Knepley      to figure out what to save where. */
115678569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
11579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1158792fecdfSBarry Smith   PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
115978569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
116078569bb4SMatthew G. Knepley     ex_block block;
116178569bb4SMatthew G. Knepley 
116278569bb4SMatthew G. Knepley     block.id   = csID[set];
116378569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1164792fecdfSBarry Smith     PetscCallExternal(ex_get_block_param, exoid, &block);
116578569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
116678569bb4SMatthew G. Knepley   }
11679566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
11689566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
11699566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
117078569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
117178569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
117278569bb4SMatthew G. Knepley 
117378569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
117478569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
117578569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
117678569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
117778569bb4SMatthew G. Knepley     if (bs == 1) {
11789566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vNatural, &varray));
1179792fecdfSBarry Smith       PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs)]);
11809566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vNatural, &varray));
118178569bb4SMatthew G. Knepley     } else {
118278569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
11839566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
11849566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
11859566063dSJacob Faibussowitsch         PetscCall(VecGetArray(vComp, &varray));
1186792fecdfSBarry Smith         PetscCallExternal(ex_get_partial_var, exoid, step, EX_ELEM_BLOCK, offset + c, csID[set], PetscMax(xs / bs - csxs, 0) + 1, csLocalSize, &varray[PetscMax(0, csxs - xs / bs)]);
11879566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(vComp, &varray));
11889566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
118978569bb4SMatthew G. Knepley       }
119078569bb4SMatthew G. Knepley     }
119178569bb4SMatthew G. Knepley     csxs += csSize[set];
119278569bb4SMatthew G. Knepley   }
11939566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
11949566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
119578569bb4SMatthew G. Knepley   if (useNatural) {
11969566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
11979566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
1198b7352c5cSAlexis Marboeuf     PetscCall(VecDestroy(&vNatural));
119978569bb4SMatthew G. Knepley   }
120078569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
120178569bb4SMatthew G. Knepley }
1202b53c8456SSatish Balay #endif
120378569bb4SMatthew G. Knepley 
12041e50132fSMatthew G. Knepley /*@
12051e50132fSMatthew G. Knepley   PetscViewerExodusIIGetId - Get the file id of the ExodusII file
12061e50132fSMatthew G. Knepley 
12071e50132fSMatthew G. Knepley   Logically Collective on PetscViewer
12081e50132fSMatthew G. Knepley 
12091e50132fSMatthew G. Knepley   Input Parameter:
12101e50132fSMatthew G. Knepley .  viewer - the PetscViewer
12111e50132fSMatthew G. Knepley 
12121e50132fSMatthew G. Knepley   Output Parameter:
1213d8d19677SJose E. Roman .  exoid - The ExodusII file id
12141e50132fSMatthew G. Knepley 
12151e50132fSMatthew G. Knepley   Level: intermediate
12161e50132fSMatthew G. Knepley 
1217db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
12181e50132fSMatthew G. Knepley @*/
12199371c9d4SSatish Balay PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid) {
12201e50132fSMatthew G. Knepley   PetscFunctionBegin;
12211e50132fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1222cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetId_C", (PetscViewer, int *), (viewer, exoid));
12236823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12246823f3c5SBlaise Bourdin }
12256823f3c5SBlaise Bourdin 
12266823f3c5SBlaise Bourdin /*@
12276823f3c5SBlaise Bourdin    PetscViewerExodusIISetOrder - Set the elements order in the exodusII file.
12286823f3c5SBlaise Bourdin 
12296823f3c5SBlaise Bourdin    Collective
12306823f3c5SBlaise Bourdin 
12316823f3c5SBlaise Bourdin    Input Parameters:
12326823f3c5SBlaise Bourdin +  viewer - the viewer
123398a6a78aSPierre Jolivet -  order - elements order
12346823f3c5SBlaise Bourdin 
12356823f3c5SBlaise Bourdin    Output Parameter:
12366823f3c5SBlaise Bourdin 
12376823f3c5SBlaise Bourdin    Level: beginner
12386823f3c5SBlaise Bourdin 
12396823f3c5SBlaise Bourdin    Note:
12406823f3c5SBlaise Bourdin 
1241db781477SPatrick Sanan .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
12426823f3c5SBlaise Bourdin @*/
12439371c9d4SSatish Balay PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order) {
12446823f3c5SBlaise Bourdin   PetscFunctionBegin;
12456823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1246cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerSetOrder_C", (PetscViewer, PetscInt), (viewer, order));
12476823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12486823f3c5SBlaise Bourdin }
12496823f3c5SBlaise Bourdin 
12506823f3c5SBlaise Bourdin /*@
12516823f3c5SBlaise Bourdin    PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file.
12526823f3c5SBlaise Bourdin 
12536823f3c5SBlaise Bourdin    Collective
12546823f3c5SBlaise Bourdin 
12556823f3c5SBlaise Bourdin    Input Parameters:
12566823f3c5SBlaise Bourdin +  viewer - the viewer
125798a6a78aSPierre Jolivet -  order - elements order
12586823f3c5SBlaise Bourdin 
12596823f3c5SBlaise Bourdin    Output Parameter:
12606823f3c5SBlaise Bourdin 
12616823f3c5SBlaise Bourdin    Level: beginner
12626823f3c5SBlaise Bourdin 
12636823f3c5SBlaise Bourdin    Note:
12646823f3c5SBlaise Bourdin 
1265db781477SPatrick Sanan .seealso: `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
12666823f3c5SBlaise Bourdin @*/
12679371c9d4SSatish Balay PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order) {
12686823f3c5SBlaise Bourdin   PetscFunctionBegin;
12696823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1270cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetOrder_C", (PetscViewer, PetscInt *), (viewer, order));
12716823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12726823f3c5SBlaise Bourdin }
12736823f3c5SBlaise Bourdin 
12746823f3c5SBlaise Bourdin /*@C
12756823f3c5SBlaise Bourdin    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
12766823f3c5SBlaise Bourdin 
12776823f3c5SBlaise Bourdin    Collective
12786823f3c5SBlaise Bourdin 
12796823f3c5SBlaise Bourdin    Input Parameters:
12806823f3c5SBlaise Bourdin +  comm - MPI communicator
12816823f3c5SBlaise Bourdin .  name - name of file
12826823f3c5SBlaise Bourdin -  type - type of file
12836823f3c5SBlaise Bourdin $    FILE_MODE_WRITE - create new file for binary output
12846823f3c5SBlaise Bourdin $    FILE_MODE_READ - open existing file for binary input
12856823f3c5SBlaise Bourdin $    FILE_MODE_APPEND - open existing file for binary output
12866823f3c5SBlaise Bourdin 
12876823f3c5SBlaise Bourdin    Output Parameter:
12886823f3c5SBlaise Bourdin .  exo - PetscViewer for Exodus II input/output to use with the specified file
12896823f3c5SBlaise Bourdin 
12906823f3c5SBlaise Bourdin    Level: beginner
12916823f3c5SBlaise Bourdin 
12926823f3c5SBlaise Bourdin    Note:
12936823f3c5SBlaise Bourdin    This PetscViewer should be destroyed with PetscViewerDestroy().
12946823f3c5SBlaise Bourdin 
1295db781477SPatrick Sanan .seealso: `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1296db781477SPatrick Sanan           `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
12976823f3c5SBlaise Bourdin @*/
12989371c9d4SSatish Balay PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo) {
12996823f3c5SBlaise Bourdin   PetscFunctionBegin;
13009566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, exo));
13019566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII));
13029566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*exo, type));
13039566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*exo, name));
13049566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*exo));
13051e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
13061e50132fSMatthew G. Knepley }
13071e50132fSMatthew G. Knepley 
130833751fbdSMichael Lange /*@C
1309eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
131033751fbdSMichael Lange 
1311d083f849SBarry Smith   Collective
131233751fbdSMichael Lange 
131333751fbdSMichael Lange   Input Parameters:
131433751fbdSMichael Lange + comm  - The MPI communicator
131533751fbdSMichael Lange . filename - The name of the ExodusII file
131633751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
131733751fbdSMichael Lange 
131833751fbdSMichael Lange   Output Parameter:
131933751fbdSMichael Lange . dm  - The DM object representing the mesh
132033751fbdSMichael Lange 
132133751fbdSMichael Lange   Level: beginner
132233751fbdSMichael Lange 
1323db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()`
132433751fbdSMichael Lange @*/
13259371c9d4SSatish Balay PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm) {
132633751fbdSMichael Lange   PetscMPIInt rank;
132733751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1328e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
132933751fbdSMichael Lange   float version;
133033751fbdSMichael Lange #endif
133133751fbdSMichael Lange 
133233751fbdSMichael Lange   PetscFunctionBegin;
133333751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
13349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
133533751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1336dd400576SPatrick Sanan   if (rank == 0) {
133733751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
133808401ef6SPierre Jolivet     PetscCheck(exoid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
133933751fbdSMichael Lange   }
13409566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm));
134148a46eb9SPierre Jolivet   if (rank == 0) PetscCallExternal(ex_close, exoid);
1342b458e8f1SJose E. Roman   PetscFunctionReturn(0);
134333751fbdSMichael Lange #else
134433751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
134533751fbdSMichael Lange #endif
134633751fbdSMichael Lange }
134733751fbdSMichael Lange 
13488f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
13499371c9d4SSatish Balay static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct) {
13508f861fbcSMatthew G. Knepley   PetscBool flg;
13518f861fbcSMatthew G. Knepley 
13528f861fbcSMatthew G. Knepley   PetscFunctionBegin;
13538f861fbcSMatthew G. Knepley   *ct = DM_POLYTOPE_UNKNOWN;
13549566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI", &flg));
13559371c9d4SSatish Balay   if (flg) {
13569371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRIANGLE;
13579371c9d4SSatish Balay     goto done;
13589371c9d4SSatish Balay   }
13599566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI3", &flg));
13609371c9d4SSatish Balay   if (flg) {
13619371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRIANGLE;
13629371c9d4SSatish Balay     goto done;
13639371c9d4SSatish Balay   }
13649566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD", &flg));
13659371c9d4SSatish Balay   if (flg) {
13669371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
13679371c9d4SSatish Balay     goto done;
13689371c9d4SSatish Balay   }
13699566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg));
13709371c9d4SSatish Balay   if (flg) {
13719371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
13729371c9d4SSatish Balay     goto done;
13739371c9d4SSatish Balay   }
13749566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg));
13759371c9d4SSatish Balay   if (flg) {
13769371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
13779371c9d4SSatish Balay     goto done;
13789371c9d4SSatish Balay   }
13799566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TETRA", &flg));
13809371c9d4SSatish Balay   if (flg) {
13819371c9d4SSatish Balay     *ct = DM_POLYTOPE_TETRAHEDRON;
13829371c9d4SSatish Balay     goto done;
13839371c9d4SSatish Balay   }
13849566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TET4", &flg));
13859371c9d4SSatish Balay   if (flg) {
13869371c9d4SSatish Balay     *ct = DM_POLYTOPE_TETRAHEDRON;
13879371c9d4SSatish Balay     goto done;
13889371c9d4SSatish Balay   }
13899566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg));
13909371c9d4SSatish Balay   if (flg) {
13919371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRI_PRISM;
13929371c9d4SSatish Balay     goto done;
13939371c9d4SSatish Balay   }
13949566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX", &flg));
13959371c9d4SSatish Balay   if (flg) {
13969371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
13979371c9d4SSatish Balay     goto done;
13989371c9d4SSatish Balay   }
13999566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX8", &flg));
14009371c9d4SSatish Balay   if (flg) {
14019371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
14029371c9d4SSatish Balay     goto done;
14039371c9d4SSatish Balay   }
14049566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg));
14059371c9d4SSatish Balay   if (flg) {
14069371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
14079371c9d4SSatish Balay     goto done;
14089371c9d4SSatish Balay   }
140998921bdaSJacob Faibussowitsch   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
14108f861fbcSMatthew G. Knepley done:
14118f861fbcSMatthew G. Knepley   PetscFunctionReturn(0);
14128f861fbcSMatthew G. Knepley }
14138f861fbcSMatthew G. Knepley #endif
14148f861fbcSMatthew G. Knepley 
1415552f7358SJed Brown /*@
141633751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
1417552f7358SJed Brown 
1418d083f849SBarry Smith   Collective
1419552f7358SJed Brown 
1420552f7358SJed Brown   Input Parameters:
1421552f7358SJed Brown + comm  - The MPI communicator
1422552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1423552f7358SJed Brown - interpolate - Create faces and edges in the mesh
1424552f7358SJed Brown 
1425552f7358SJed Brown   Output Parameter:
1426552f7358SJed Brown . dm  - The DM object representing the mesh
1427552f7358SJed Brown 
1428552f7358SJed Brown   Level: beginner
1429552f7358SJed Brown 
1430db781477SPatrick Sanan .seealso: `DMPLEX`, `DMCreate()`
1431552f7358SJed Brown @*/
14329371c9d4SSatish Balay PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm) {
1433552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1434552f7358SJed Brown   PetscMPIInt  num_proc, rank;
1435ae9eebabSLisandro Dalcin   DMLabel      cellSets = NULL, faceSets = NULL, vertSets = NULL;
1436552f7358SJed Brown   PetscSection coordSection;
1437552f7358SJed Brown   Vec          coordinates;
1438552f7358SJed Brown   PetscScalar *coords;
1439552f7358SJed Brown   PetscInt     coordSize, v;
1440552f7358SJed Brown   /* Read from ex_get_init() */
1441552f7358SJed Brown   char         title[PETSC_MAX_PATH_LEN + 1];
14425f80ce2aSJacob Faibussowitsch   int          dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0;
1443552f7358SJed Brown   int          num_cs = 0, num_vs = 0, num_fs = 0;
1444552f7358SJed Brown #endif
1445552f7358SJed Brown 
1446552f7358SJed Brown   PetscFunctionBegin;
1447552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
14489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
14499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
14509566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
14519566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1452a5b23f4aSJose E. Roman   /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */
1453dd400576SPatrick Sanan   if (rank == 0) {
14549566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1));
1455792fecdfSBarry Smith     PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
145628b400f6SJacob Faibussowitsch     PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set");
1457552f7358SJed Brown   }
14589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm));
14599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
14609566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, title));
14619566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
1462412e9a14SMatthew G. Knepley   /*   We do not want this label automatically computed, instead we compute it here */
14639566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "celltype"));
1464552f7358SJed Brown 
1465552f7358SJed Brown   /* Read cell sets information */
1466dd400576SPatrick Sanan   if (rank == 0) {
1467552f7358SJed Brown     PetscInt *cone;
1468e8f6893fSMatthew G. Knepley     int       c, cs, ncs, c_loc, v, v_loc;
1469552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
1470e8f6893fSMatthew G. Knepley     int      *cs_id, *cs_order;
1471552f7358SJed Brown     /* Read from ex_get_elem_block() */
1472552f7358SJed Brown     char      buffer[PETSC_MAX_PATH_LEN + 1];
1473e8f6893fSMatthew G. Knepley     int       num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1474552f7358SJed Brown     /* Read from ex_get_elem_conn() */
1475552f7358SJed Brown     int      *cs_connect;
1476552f7358SJed Brown 
1477552f7358SJed Brown     /* Get cell sets IDs */
14789566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order));
1479792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id);
1480552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
1481552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1482e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
1483e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
14848f861fbcSMatthew G. Knepley       DMPolytopeType ct;
14858f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
14868f861fbcSMatthew G. Knepley 
14879566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1488792fecdfSBarry Smith       PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
14899566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
14908f861fbcSMatthew G. Knepley       dim = PetscMax(dim, DMPolytopeTypeGetDim(ct));
1491792fecdfSBarry Smith       PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
14928f861fbcSMatthew G. Knepley       switch (ct) {
14938f861fbcSMatthew G. Knepley       case DM_POLYTOPE_TRI_PRISM:
1494e8f6893fSMatthew G. Knepley         cs_order[cs] = cs;
1495e8f6893fSMatthew G. Knepley         ++num_hybrid;
1496e8f6893fSMatthew G. Knepley         break;
1497e8f6893fSMatthew G. Knepley       default:
1498e8f6893fSMatthew G. Knepley         for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1];
1499e8f6893fSMatthew G. Knepley         cs_order[cs - num_hybrid] = cs;
1500e8f6893fSMatthew G. Knepley       }
1501e8f6893fSMatthew G. Knepley     }
1502552f7358SJed Brown     /* First set sizes */
1503e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
15048f861fbcSMatthew G. Knepley       DMPolytopeType ct;
15058f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
1506e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
15078f861fbcSMatthew G. Knepley 
15089566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1509792fecdfSBarry Smith       PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
15109566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
1511792fecdfSBarry Smith       PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
1512552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
15139566063dSJacob Faibussowitsch         PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell));
15149566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCellType(*dm, c, ct));
1515552f7358SJed Brown       }
1516552f7358SJed Brown     }
15179566063dSJacob Faibussowitsch     for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
15189566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dm));
1519e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1520e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
1521792fecdfSBarry Smith       PetscCallExternal(ex_get_block, exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
15229566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone));
1523792fecdfSBarry Smith       PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL);
1524eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1525552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1526636e64ffSStefano Zampini         DMPolytopeType ct;
1527636e64ffSStefano Zampini 
1528ad540459SPierre Jolivet         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1;
15299566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(*dm, c, &ct));
15309566063dSJacob Faibussowitsch         PetscCall(DMPlexInvertCell(ct, cone));
15319566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCone(*dm, c, cone));
15329566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]));
1533552f7358SJed Brown       }
15349566063dSJacob Faibussowitsch       PetscCall(PetscFree2(cs_connect, cone));
1535552f7358SJed Brown     }
15369566063dSJacob Faibussowitsch     PetscCall(PetscFree2(cs_id, cs_order));
1537552f7358SJed Brown   }
15388f861fbcSMatthew G. Knepley   {
15398f861fbcSMatthew G. Knepley     PetscInt ints[] = {dim, dimEmbed};
15408f861fbcSMatthew G. Knepley 
15419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm));
15429566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dm, ints[0]));
15439566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(*dm, ints[1]));
15448f861fbcSMatthew G. Knepley     dim      = ints[0];
15458f861fbcSMatthew G. Knepley     dimEmbed = ints[1];
15468f861fbcSMatthew G. Knepley   }
15479566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
15489566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
1549552f7358SJed Brown   if (interpolate) {
15505fd9971aSMatthew G. Knepley     DM idm;
1551552f7358SJed Brown 
15529566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
15539566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1554552f7358SJed Brown     *dm = idm;
1555552f7358SJed Brown   }
1556552f7358SJed Brown 
1557552f7358SJed Brown   /* Create vertex set label */
1558dd400576SPatrick Sanan   if (rank == 0 && (num_vs > 0)) {
1559552f7358SJed Brown     int  vs, v;
1560552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1561552f7358SJed Brown     int *vs_id;
1562552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1563f958083aSBlaise Bourdin     int  num_vertex_in_set;
1564552f7358SJed Brown     /* Read from ex_get_node_set() */
1565552f7358SJed Brown     int *vs_vertex_list;
1566552f7358SJed Brown 
1567552f7358SJed Brown     /* Get vertex set ids */
15689566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_vs, &vs_id));
1569792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id);
1570552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
1571792fecdfSBarry Smith       PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
15729566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list));
1573792fecdfSBarry Smith       PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);
157448a46eb9SPierre Jolivet       for (v = 0; v < num_vertex_in_set; ++v) PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v] + numCells - 1, vs_id[vs]));
15759566063dSJacob Faibussowitsch       PetscCall(PetscFree(vs_vertex_list));
1576552f7358SJed Brown     }
15779566063dSJacob Faibussowitsch     PetscCall(PetscFree(vs_id));
1578552f7358SJed Brown   }
1579552f7358SJed Brown   /* Read coordinates */
15809566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
15819566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
15829566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
15839566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
1584552f7358SJed Brown   for (v = numCells; v < numCells + numVertices; ++v) {
15859566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
15869566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
1587552f7358SJed Brown   }
15889566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
15899566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
15909566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
15919566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
15929566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
15939566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, dimEmbed));
15949566063dSJacob Faibussowitsch   PetscCall(VecSetType(coordinates, VECSTANDARD));
15959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
1596dd400576SPatrick Sanan   if (rank == 0) {
1597e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1598552f7358SJed Brown 
15999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z));
1600792fecdfSBarry Smith     PetscCallExternal(ex_get_coord, exoid, x, y, z);
16018f861fbcSMatthew G. Knepley     if (dimEmbed > 0) {
16028f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v];
16030d644c17SKarl Rupp     }
16048f861fbcSMatthew G. Knepley     if (dimEmbed > 1) {
16058f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v];
16060d644c17SKarl Rupp     }
16078f861fbcSMatthew G. Knepley     if (dimEmbed > 2) {
16088f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v];
16090d644c17SKarl Rupp     }
16109566063dSJacob Faibussowitsch     PetscCall(PetscFree3(x, y, z));
1611552f7358SJed Brown   }
16129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
16139566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
16149566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
1615552f7358SJed Brown 
1616552f7358SJed Brown   /* Create side set label */
1617dd400576SPatrick Sanan   if (rank == 0 && interpolate && (num_fs > 0)) {
1618552f7358SJed Brown     int    fs, f, voff;
1619552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1620552f7358SJed Brown     int   *fs_id;
1621552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1622f958083aSBlaise Bourdin     int    num_side_in_set;
1623552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1624552f7358SJed Brown     int   *fs_vertex_count_list, *fs_vertex_list;
1625ef073283Smichael_afanasiev     /* Read side set labels */
1626ef073283Smichael_afanasiev     char   fs_name[MAX_STR_LENGTH + 1];
1627034d25ccSMatthew G. Knepley     size_t fs_name_len;
1628552f7358SJed Brown 
1629552f7358SJed Brown     /* Get side set ids */
16309566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_fs, &fs_id));
1631792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id);
1632552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
1633792fecdfSBarry Smith       PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);
16349566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list));
1635792fecdfSBarry Smith       PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);
1636ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1637ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1638034d25ccSMatthew G. Knepley       if (!fs_name_err) {
1639034d25ccSMatthew G. Knepley         PetscCall(PetscStrlen(fs_name, &fs_name_len));
1640034d25ccSMatthew G. Knepley         if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH));
1641034d25ccSMatthew G. Knepley       }
1642552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
16430298fd71SBarry Smith         const PetscInt *faces    = NULL;
1644552f7358SJed Brown         PetscInt        faceSize = fs_vertex_count_list[f], numFaces;
1645552f7358SJed Brown         PetscInt        faceVertices[4], v;
1646552f7358SJed Brown 
1647197f6770SJed Brown         PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize);
1648ad540459SPierre Jolivet         for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1;
16499566063dSJacob Faibussowitsch         PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1650197f6770SJed Brown         PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces);
16519566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]));
1652ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1653034d25ccSMatthew G. Knepley         if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]));
16549566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1655552f7358SJed Brown       }
16569566063dSJacob Faibussowitsch       PetscCall(PetscFree2(fs_vertex_count_list, fs_vertex_list));
1657552f7358SJed Brown     }
16589566063dSJacob Faibussowitsch     PetscCall(PetscFree(fs_id));
1659552f7358SJed Brown   }
1660ae9eebabSLisandro Dalcin 
1661ae9eebabSLisandro Dalcin   { /* Create Cell/Face/Vertex Sets labels at all processes */
16629371c9d4SSatish Balay     enum {
16639371c9d4SSatish Balay       n = 3
16649371c9d4SSatish Balay     };
1665ae9eebabSLisandro Dalcin     PetscBool flag[n];
1666ae9eebabSLisandro Dalcin 
1667ae9eebabSLisandro Dalcin     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1668ae9eebabSLisandro Dalcin     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1669ae9eebabSLisandro Dalcin     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
16709566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
16719566063dSJacob Faibussowitsch     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
16729566063dSJacob Faibussowitsch     if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
16739566063dSJacob Faibussowitsch     if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
1674ae9eebabSLisandro Dalcin   }
1675b458e8f1SJose E. Roman   PetscFunctionReturn(0);
1676552f7358SJed Brown #else
1677552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1678552f7358SJed Brown #endif
1679552f7358SJed Brown }
1680