xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision bdbfaa5d8000d6881a6854c721cdde3d6a8791e9)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown 
4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII)
5552f7358SJed Brown   #include <netcdf.h>
6552f7358SJed Brown   #include <exodusII.h>
7552f7358SJed Brown #endif
8552f7358SJed Brown 
91e50132fSMatthew G. Knepley #include <petsc/private/viewerimpl.h>
101e50132fSMatthew G. Knepley #include <petsc/private/viewerexodusiiimpl.h>
11e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
12a1cb98faSBarry Smith /*@C
13a1cb98faSBarry Smith   PETSC_VIEWER_EXODUSII_ - Creates an `PETSCVIEWEREXODUSII` `PetscViewer` shared by all processors in a communicator.
141e50132fSMatthew G. Knepley 
151e50132fSMatthew G. Knepley   Collective
161e50132fSMatthew G. Knepley 
171e50132fSMatthew G. Knepley   Input Parameter:
18a1cb98faSBarry Smith . comm - the MPI communicator to share the `PETSCVIEWEREXODUSII` `PetscViewer`
191e50132fSMatthew G. Knepley 
201e50132fSMatthew G. Knepley   Level: intermediate
211e50132fSMatthew G. Knepley 
22a1cb98faSBarry Smith   Note:
231e50132fSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return
241e50132fSMatthew G. Knepley   an error code.  The GLVIS PetscViewer is usually used in the form
251e50132fSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
261e50132fSMatthew G. Knepley 
27a1cb98faSBarry Smith   Fortran Note:
28a1cb98faSBarry Smith   No support in Fortran
29a1cb98faSBarry Smith 
30a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewer`, `PetscViewerExodusIIOpen()`, `PetscViewerType`, `PetscViewerCreate()`, `PetscViewerDestroy()`
31a1cb98faSBarry Smith @*/
32d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm)
33d71ae5a4SJacob Faibussowitsch {
341e50132fSMatthew G. Knepley   PetscViewer    viewer;
351e50132fSMatthew G. Knepley   PetscErrorCode ierr;
361e50132fSMatthew G. Knepley 
371e50132fSMatthew G. Knepley   PetscFunctionBegin;
381e50132fSMatthew G. Knepley   ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer);
399371c9d4SSatish Balay   if (ierr) {
409371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
4111cc89d2SBarry Smith     PetscFunctionReturn(NULL);
429371c9d4SSatish Balay   }
431e50132fSMatthew G. Knepley   ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
449371c9d4SSatish Balay   if (ierr) {
459371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_EXODUSII_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
4611cc89d2SBarry Smith     PetscFunctionReturn(NULL);
479371c9d4SSatish Balay   }
481e50132fSMatthew G. Knepley   PetscFunctionReturn(viewer);
491e50132fSMatthew G. Knepley }
501e50132fSMatthew G. Knepley 
51d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer)
52d71ae5a4SJacob Faibussowitsch {
536823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)v->data;
541e50132fSMatthew G. Knepley 
551e50132fSMatthew G. Knepley   PetscFunctionBegin;
569566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename:    %s\n", exo->filename));
579566063dSJacob Faibussowitsch   if (exo->exoid) PetscCall(PetscViewerASCIIPrintf(viewer, "exoid:       %d\n", exo->exoid));
589566063dSJacob Faibussowitsch   if (exo->btype) PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode:     %d\n", exo->btype));
59197f6770SJed Brown   if (exo->order) PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order:  %" PetscInt_FMT "\n", exo->order));
601e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
611e50132fSMatthew G. Knepley }
621e50132fSMatthew G. Knepley 
63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscViewer v, PetscOptionItems *PetscOptionsObject)
64d71ae5a4SJacob Faibussowitsch {
651e50132fSMatthew G. Knepley   PetscFunctionBegin;
66d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "ExodusII PetscViewer Options");
67d0609cedSBarry Smith   PetscOptionsHeadEnd();
681e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
691e50132fSMatthew G. Knepley }
701e50132fSMatthew G. Knepley 
71d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer)
72d71ae5a4SJacob Faibussowitsch {
731e50132fSMatthew G. Knepley   PetscFunctionBegin;
741e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
751e50132fSMatthew G. Knepley }
761e50132fSMatthew G. Knepley 
77d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
78d71ae5a4SJacob Faibussowitsch {
791e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
801e50132fSMatthew G. Knepley 
811e50132fSMatthew G. Knepley   PetscFunctionBegin;
8248a46eb9SPierre Jolivet   if (exo->exoid >= 0) PetscCallExternal(ex_close, exo->exoid);
839566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo->filename));
849566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo));
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
882e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetId_C", NULL));
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerGetOrder_C", NULL));
919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerSetOrder_C", NULL));
921e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
931e50132fSMatthew G. Knepley }
941e50132fSMatthew G. Knepley 
95d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
96d71ae5a4SJacob Faibussowitsch {
971e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
981e50132fSMatthew G. Knepley   PetscMPIInt           rank;
991e50132fSMatthew G. Knepley   int                   CPU_word_size, IO_word_size, EXO_mode;
1006823f3c5SBlaise Bourdin   MPI_Info              mpi_info = MPI_INFO_NULL;
1016823f3c5SBlaise Bourdin   float                 EXO_version;
1021e50132fSMatthew G. Knepley 
1039566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
1041e50132fSMatthew G. Knepley   CPU_word_size = sizeof(PetscReal);
1051e50132fSMatthew G. Knepley   IO_word_size  = sizeof(PetscReal);
1061e50132fSMatthew G. Knepley 
1071e50132fSMatthew G. Knepley   PetscFunctionBegin;
1086823f3c5SBlaise Bourdin   if (exo->exoid >= 0) {
109792fecdfSBarry Smith     PetscCallExternal(ex_close, exo->exoid);
1106823f3c5SBlaise Bourdin     exo->exoid = -1;
1116823f3c5SBlaise Bourdin   }
1129566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscFree(exo->filename));
1139566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &exo->filename));
1141e50132fSMatthew G. Knepley   switch (exo->btype) {
115d71ae5a4SJacob Faibussowitsch   case FILE_MODE_READ:
116d71ae5a4SJacob Faibussowitsch     EXO_mode = EX_READ;
117d71ae5a4SJacob Faibussowitsch     break;
1181e50132fSMatthew G. Knepley   case FILE_MODE_APPEND:
1196823f3c5SBlaise Bourdin   case FILE_MODE_UPDATE:
1206823f3c5SBlaise Bourdin   case FILE_MODE_APPEND_UPDATE:
1216823f3c5SBlaise Bourdin     /* Will fail if the file does not already exist */
1226823f3c5SBlaise Bourdin     EXO_mode = EX_WRITE;
1231e50132fSMatthew G. Knepley     break;
1241e50132fSMatthew G. Knepley   case FILE_MODE_WRITE:
1256823f3c5SBlaise Bourdin     /*
1266823f3c5SBlaise Bourdin       exodus only allows writing geometry upon file creation, so we will let DMView create the file.
1276823f3c5SBlaise Bourdin     */
1286823f3c5SBlaise Bourdin     PetscFunctionReturn(0);
1291e50132fSMatthew G. Knepley     break;
130d71ae5a4SJacob Faibussowitsch   default:
131d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
1321e50132fSMatthew G. Knepley   }
1331e50132fSMatthew G. Knepley   #if defined(PETSC_USE_64BIT_INDICES)
1341e50132fSMatthew G. Knepley   EXO_mode += EX_ALL_INT64_API;
1351e50132fSMatthew G. Knepley   #endif
1366823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(name, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, PETSC_COMM_WORLD, mpi_info);
13708401ef6SPierre Jolivet   PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
1381e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1391e50132fSMatthew G. Knepley }
1401e50132fSMatthew G. Knepley 
141d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name)
142d71ae5a4SJacob Faibussowitsch {
1431e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1441e50132fSMatthew G. Knepley 
1451e50132fSMatthew G. Knepley   PetscFunctionBegin;
1461e50132fSMatthew G. Knepley   *name = exo->filename;
1471e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1481e50132fSMatthew G. Knepley }
1491e50132fSMatthew G. Knepley 
150d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
151d71ae5a4SJacob Faibussowitsch {
1521e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1531e50132fSMatthew G. Knepley 
1541e50132fSMatthew G. Knepley   PetscFunctionBegin;
1551e50132fSMatthew G. Knepley   exo->btype = type;
1561e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1571e50132fSMatthew G. Knepley }
1581e50132fSMatthew G. Knepley 
159d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
160d71ae5a4SJacob Faibussowitsch {
1611e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1621e50132fSMatthew G. Knepley 
1631e50132fSMatthew G. Knepley   PetscFunctionBegin;
1641e50132fSMatthew G. Knepley   *type = exo->btype;
1651e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1661e50132fSMatthew G. Knepley }
1671e50132fSMatthew G. Knepley 
168d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid)
169d71ae5a4SJacob Faibussowitsch {
1701e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1711e50132fSMatthew G. Knepley 
1721e50132fSMatthew G. Knepley   PetscFunctionBegin;
1731e50132fSMatthew G. Knepley   *exoid = exo->exoid;
1741e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1751e50132fSMatthew G. Knepley }
1761e50132fSMatthew G. Knepley 
177d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order)
178d71ae5a4SJacob Faibussowitsch {
1796823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1801e50132fSMatthew G. Knepley 
1811e50132fSMatthew G. Knepley   PetscFunctionBegin;
1826823f3c5SBlaise Bourdin   *order = exo->order;
1836823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
1846823f3c5SBlaise Bourdin }
1856823f3c5SBlaise Bourdin 
186d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order)
187d71ae5a4SJacob Faibussowitsch {
1886823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
1896823f3c5SBlaise Bourdin 
1906823f3c5SBlaise Bourdin   PetscFunctionBegin;
1916823f3c5SBlaise Bourdin   exo->order = order;
1921e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1931e50132fSMatthew G. Knepley }
1941e50132fSMatthew G. Knepley 
1951e50132fSMatthew G. Knepley /*MC
19600373969SVaclav Hapla    PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
1971e50132fSMatthew G. Knepley 
198db781477SPatrick Sanan .seealso: `PetscViewerExodusIIOpen()`, `PetscViewerCreate()`, `PETSCVIEWERBINARY`, `PETSCVIEWERHDF5`, `DMView()`,
199db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
2001e50132fSMatthew G. Knepley 
2011e50132fSMatthew G. Knepley   Level: beginner
2021e50132fSMatthew G. Knepley M*/
2031e50132fSMatthew G. Knepley 
204d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
205d71ae5a4SJacob Faibussowitsch {
2061e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo;
2071e50132fSMatthew G. Knepley 
2081e50132fSMatthew G. Knepley   PetscFunctionBegin;
2094dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&exo));
2101e50132fSMatthew G. Knepley 
2111e50132fSMatthew G. Knepley   v->data                = (void *)exo;
2121e50132fSMatthew G. Knepley   v->ops->destroy        = PetscViewerDestroy_ExodusII;
2131e50132fSMatthew G. Knepley   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
2141e50132fSMatthew G. Knepley   v->ops->setup          = PetscViewerSetUp_ExodusII;
2151e50132fSMatthew G. Knepley   v->ops->view           = PetscViewerView_ExodusII;
2161e50132fSMatthew G. Knepley   v->ops->flush          = 0;
2171e50132fSMatthew G. Knepley   exo->btype             = (PetscFileMode)-1;
2181e50132fSMatthew G. Knepley   exo->filename          = 0;
2191e50132fSMatthew G. Knepley   exo->exoid             = -1;
2201e50132fSMatthew G. Knepley 
2219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_ExodusII));
2229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_ExodusII));
2239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_ExodusII));
2249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_ExodusII));
2259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetId_C", PetscViewerExodusIIGetId_ExodusII));
2269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerSetOrder_C", PetscViewerExodusIISetOrder_ExodusII));
2279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerGetOrder_C", PetscViewerExodusIIGetOrder_ExodusII));
2281e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
2291e50132fSMatthew G. Knepley }
2301e50132fSMatthew G. Knepley 
2311e50132fSMatthew G. Knepley /*
23278569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
23378569bb4SMatthew G. Knepley 
2346823f3c5SBlaise Bourdin   Collective
23578569bb4SMatthew G. Knepley 
23678569bb4SMatthew G. Knepley   Input Parameters:
23778569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
23878569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
23978569bb4SMatthew G. Knepley - name     - the name of the result
24078569bb4SMatthew G. Knepley 
24178569bb4SMatthew G. Knepley   Output Parameters:
24278569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
24378569bb4SMatthew G. Knepley 
24478569bb4SMatthew G. Knepley   Notes:
24578569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
24678569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
24778569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
24878569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
24978569bb4SMatthew G. Knepley 
25078569bb4SMatthew G. Knepley   Level: beginner
25178569bb4SMatthew G. Knepley 
252c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadNodal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
25378569bb4SMatthew G. Knepley */
254d71ae5a4SJacob Faibussowitsch PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
255d71ae5a4SJacob Faibussowitsch {
2566de834b4SMatthew G. Knepley   int       num_vars, i, j;
25778569bb4SMatthew G. Knepley   char      ext_name[MAX_STR_LENGTH + 1], var_name[MAX_STR_LENGTH + 1];
25878569bb4SMatthew G. Knepley   const int num_suffix = 5;
25978569bb4SMatthew G. Knepley   char     *suffix[5];
26078569bb4SMatthew G. Knepley   PetscBool flg;
26178569bb4SMatthew G. Knepley 
26278569bb4SMatthew G. Knepley   PetscFunctionBegin;
26378569bb4SMatthew G. Knepley   suffix[0] = (char *)"";
26478569bb4SMatthew G. Knepley   suffix[1] = (char *)"_X";
26578569bb4SMatthew G. Knepley   suffix[2] = (char *)"_XX";
26678569bb4SMatthew G. Knepley   suffix[3] = (char *)"_1";
26778569bb4SMatthew G. Knepley   suffix[4] = (char *)"_11";
26878569bb4SMatthew G. Knepley 
2696823f3c5SBlaise Bourdin   *varIndex = -1;
270792fecdfSBarry Smith   PetscCallExternal(ex_get_variable_param, exoid, obj_type, &num_vars);
27178569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
272792fecdfSBarry Smith     PetscCallExternal(ex_get_variable_name, exoid, obj_type, i + 1, var_name);
27378569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j) {
2749566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
2759566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
2769566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(ext_name, var_name, &flg));
277ad540459SPierre Jolivet       if (flg) *varIndex = i + 1;
2786823f3c5SBlaise Bourdin     }
2796823f3c5SBlaise Bourdin   }
28078569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
28178569bb4SMatthew G. Knepley }
28278569bb4SMatthew G. Knepley 
28378569bb4SMatthew G. Knepley /*
2846823f3c5SBlaise Bourdin   DMView_PlexExodusII - Write a DM to disk in exodus format
28578569bb4SMatthew G. Knepley 
28678569bb4SMatthew G. Knepley   Collective on dm
28778569bb4SMatthew G. Knepley 
28878569bb4SMatthew G. Knepley   Input Parameters:
28978569bb4SMatthew G. Knepley + dm  - The dm to be written
2906823f3c5SBlaise Bourdin . viewer - an exodusII viewer
29178569bb4SMatthew G. Knepley 
29278569bb4SMatthew G. Knepley   Notes:
29378569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
29478569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
29578569bb4SMatthew G. Knepley 
2966823f3c5SBlaise Bourdin   If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices)
2976823f3c5SBlaise Bourdin   will be written.
29878569bb4SMatthew G. Knepley 
29978569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
30078569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
3016823f3c5SBlaise Bourdin   The order of the mesh shall be set using PetscViewerExodusIISetOrder
30278569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
30378569bb4SMatthew G. Knepley 
3046823f3c5SBlaise Bourdin   This function will only handle TRI, TET, QUAD, and HEX cells.
30578569bb4SMatthew G. Knepley   Level: beginner
30678569bb4SMatthew G. Knepley 
3076823f3c5SBlaise Bourdin .seealso:
30878569bb4SMatthew G. Knepley */
309d71ae5a4SJacob Faibussowitsch PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer)
310d71ae5a4SJacob Faibussowitsch {
3119371c9d4SSatish Balay   enum ElemType {
3129371c9d4SSatish Balay     TRI,
3139371c9d4SSatish Balay     QUAD,
3149371c9d4SSatish Balay     TET,
3159371c9d4SSatish Balay     HEX
3169371c9d4SSatish Balay   };
31778569bb4SMatthew G. Knepley   MPI_Comm comm;
3186823f3c5SBlaise Bourdin   PetscInt degree; /* the order of the mesh */
31978569bb4SMatthew G. Knepley   /* Connectivity Variables */
32078569bb4SMatthew G. Knepley   PetscInt cellsNotInConnectivity;
32178569bb4SMatthew G. Knepley   /* Cell Sets */
32278569bb4SMatthew G. Knepley   DMLabel         csLabel;
32378569bb4SMatthew G. Knepley   IS              csIS;
32478569bb4SMatthew G. Knepley   const PetscInt *csIdx;
32578569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
32678569bb4SMatthew G. Knepley   enum ElemType  *type;
327fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
32878569bb4SMatthew G. Knepley   /* Coordinate Variables */
32978569bb4SMatthew G. Knepley   DM           cdm;
3306823f3c5SBlaise Bourdin   PetscSection coordSection;
33178569bb4SMatthew G. Knepley   Vec          coord;
33278569bb4SMatthew G. Knepley   PetscInt   **nodes;
333eae8a387SMatthew G. Knepley   PetscInt     depth, d, dim, skipCells = 0;
33478569bb4SMatthew G. Knepley   PetscInt     pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
33578569bb4SMatthew G. Knepley   PetscInt     num_vs, num_fs;
33678569bb4SMatthew G. Knepley   PetscMPIInt  rank, size;
33778569bb4SMatthew G. Knepley   const char  *dmName;
338fe209ef9SBlaise Bourdin   PetscInt     nodesTriP1[4]  = {3, 0, 0, 0};
339fe209ef9SBlaise Bourdin   PetscInt     nodesTriP2[4]  = {3, 3, 0, 0};
340fe209ef9SBlaise Bourdin   PetscInt     nodesQuadP1[4] = {4, 0, 0, 0};
341fe209ef9SBlaise Bourdin   PetscInt     nodesQuadP2[4] = {4, 4, 0, 1};
342fe209ef9SBlaise Bourdin   PetscInt     nodesTetP1[4]  = {4, 0, 0, 0};
343fe209ef9SBlaise Bourdin   PetscInt     nodesTetP2[4]  = {4, 6, 0, 0};
344fe209ef9SBlaise Bourdin   PetscInt     nodesHexP1[4]  = {8, 0, 0, 0};
345fe209ef9SBlaise Bourdin   PetscInt     nodesHexP2[4]  = {8, 12, 6, 1};
3466823f3c5SBlaise Bourdin   int          CPU_word_size, IO_word_size, EXO_mode;
3476823f3c5SBlaise Bourdin   float        EXO_version;
3486823f3c5SBlaise Bourdin 
3496823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *)viewer->data;
3506823f3c5SBlaise Bourdin 
35178569bb4SMatthew G. Knepley   PetscFunctionBegin;
3529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3556823f3c5SBlaise Bourdin 
3566823f3c5SBlaise Bourdin   /*
3576823f3c5SBlaise Bourdin     Creating coordSection is a collective operation so we do it somewhat out of sequence
3586823f3c5SBlaise Bourdin   */
3599566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &coordSection));
3609566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
361*bdbfaa5dSBlaise Bourdin   /*
362*bdbfaa5dSBlaise Bourdin     Check that all points are on rank 0 since we don't know how to save distributed DM in exodus format
363*bdbfaa5dSBlaise Bourdin   */
364*bdbfaa5dSBlaise Bourdin   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
365*bdbfaa5dSBlaise Bourdin   PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
366*bdbfaa5dSBlaise Bourdin   PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
367*bdbfaa5dSBlaise Bourdin   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
368*bdbfaa5dSBlaise Bourdin   numCells    = cEnd - cStart;
369*bdbfaa5dSBlaise Bourdin   numEdges    = eEnd - eStart;
370*bdbfaa5dSBlaise Bourdin   numVertices = vEnd - vStart;
371*bdbfaa5dSBlaise Bourdin   PetscCheck(!(rank && (numCells || numEdges || numVertices)), PETSC_COMM_SELF, PETSC_ERR_SUP, "Writing distributed DM in exodusII format not supported");
372c5853193SPierre Jolivet   if (rank == 0) {
3736823f3c5SBlaise Bourdin     switch (exo->btype) {
3746823f3c5SBlaise Bourdin     case FILE_MODE_READ:
3756823f3c5SBlaise Bourdin     case FILE_MODE_APPEND:
3766823f3c5SBlaise Bourdin     case FILE_MODE_UPDATE:
3776823f3c5SBlaise Bourdin     case FILE_MODE_APPEND_UPDATE:
3786823f3c5SBlaise Bourdin       /* exodusII does not allow writing geometry to an existing file */
37998921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
3806823f3c5SBlaise Bourdin     case FILE_MODE_WRITE:
3816823f3c5SBlaise Bourdin       /* Create an empty file if one already exists*/
3826823f3c5SBlaise Bourdin       EXO_mode = EX_CLOBBER;
3836823f3c5SBlaise Bourdin   #if defined(PETSC_USE_64BIT_INDICES)
3846823f3c5SBlaise Bourdin       EXO_mode += EX_ALL_INT64_API;
3856823f3c5SBlaise Bourdin   #endif
3866823f3c5SBlaise Bourdin       CPU_word_size = sizeof(PetscReal);
3876823f3c5SBlaise Bourdin       IO_word_size  = sizeof(PetscReal);
3886823f3c5SBlaise Bourdin       exo->exoid    = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
38908401ef6SPierre Jolivet       PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
3906823f3c5SBlaise Bourdin 
3916823f3c5SBlaise Bourdin       break;
392d71ae5a4SJacob Faibussowitsch     default:
393d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3946823f3c5SBlaise Bourdin     }
3956823f3c5SBlaise Bourdin 
39678569bb4SMatthew G. Knepley     /* --- Get DM info --- */
3979566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject)dm, &dmName));
3989566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepth(dm, &depth));
3999566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
4009566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4019371c9d4SSatish Balay     if (depth == 3) {
4029371c9d4SSatish Balay       numFaces = fEnd - fStart;
4039371c9d4SSatish Balay     } else {
4049371c9d4SSatish Balay       numFaces = 0;
4059371c9d4SSatish Balay     }
4069566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs));
4079566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs));
4089566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs));
4099566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coord));
4109566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
41178569bb4SMatthew G. Knepley     if (num_cs > 0) {
4129566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel));
4139566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(csLabel, &csIS));
4149566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(csIS, &csIdx));
41578569bb4SMatthew G. Knepley     }
4169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &nodes));
41778569bb4SMatthew G. Knepley     /* Set element type for each block and compute total number of nodes */
4189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &type));
41978569bb4SMatthew G. Knepley     numNodes = numVertices;
4206823f3c5SBlaise Bourdin 
4219566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIGetOrder(viewer, &degree));
422ad540459SPierre Jolivet     if (degree == 2) numNodes += numEdges;
42378569bb4SMatthew G. Knepley     cellsNotInConnectivity = numCells;
42478569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
42578569bb4SMatthew G. Knepley       IS              stratumIS;
42678569bb4SMatthew G. Knepley       const PetscInt *cells;
42778569bb4SMatthew G. Knepley       PetscScalar    *xyz = NULL;
42878569bb4SMatthew G. Knepley       PetscInt        csSize, closureSize;
42978569bb4SMatthew G. Knepley 
4309566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4319566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4329566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
4339566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
43478569bb4SMatthew G. Knepley       switch (dim) {
43578569bb4SMatthew G. Knepley       case 2:
4369371c9d4SSatish Balay         if (closureSize == 3 * dim) {
4379371c9d4SSatish Balay           type[cs] = TRI;
4389371c9d4SSatish Balay         } else if (closureSize == 4 * dim) {
4399371c9d4SSatish Balay           type[cs] = QUAD;
4409371c9d4SSatish 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);
44178569bb4SMatthew G. Knepley         break;
44278569bb4SMatthew G. Knepley       case 3:
4439371c9d4SSatish Balay         if (closureSize == 4 * dim) {
4449371c9d4SSatish Balay           type[cs] = TET;
4459371c9d4SSatish Balay         } else if (closureSize == 8 * dim) {
4469371c9d4SSatish Balay           type[cs] = HEX;
4479371c9d4SSatish 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);
44878569bb4SMatthew G. Knepley         break;
449d71ae5a4SJacob Faibussowitsch       default:
450d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not handled by ExodusII viewer", dim);
45178569bb4SMatthew G. Knepley       }
452ad540459SPierre Jolivet       if ((degree == 2) && (type[cs] == QUAD)) numNodes += csSize;
4539371c9d4SSatish Balay       if ((degree == 2) && (type[cs] == HEX)) {
4549371c9d4SSatish Balay         numNodes += csSize;
4559371c9d4SSatish Balay         numNodes += numFaces;
4569371c9d4SSatish Balay       }
4579566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
45878569bb4SMatthew G. Knepley       /* Set nodes and Element type */
45978569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
46078569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesTriP1;
46178569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTriP2;
46278569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
46378569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesQuadP1;
46478569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesQuadP2;
46578569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
46678569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesTetP1;
46778569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTetP2;
46878569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
46978569bb4SMatthew G. Knepley         if (degree == 1) nodes[cs] = nodesHexP1;
47078569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesHexP2;
47178569bb4SMatthew G. Knepley       }
47278569bb4SMatthew G. Knepley       /* Compute the number of cells not in the connectivity table */
47378569bb4SMatthew G. Knepley       cellsNotInConnectivity -= nodes[cs][3] * csSize;
47478569bb4SMatthew G. Knepley 
4759566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
4769566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
47778569bb4SMatthew G. Knepley     }
478*bdbfaa5dSBlaise Bourdin     if (num_cs) PetscCallExternal(ex_put_init, exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);
47978569bb4SMatthew G. Knepley     /* --- Connectivity --- */
48078569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
48178569bb4SMatthew G. Knepley       IS              stratumIS;
48278569bb4SMatthew G. Knepley       const PetscInt *cells;
483fe209ef9SBlaise Bourdin       PetscInt       *connect, off = 0;
484eae8a387SMatthew G. Knepley       PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
48578569bb4SMatthew G. Knepley       PetscInt        csSize, c, connectSize, closureSize;
48678569bb4SMatthew G. Knepley       char           *elem_type        = NULL;
48778569bb4SMatthew G. Knepley       char            elem_type_tri3[] = "TRI3", elem_type_quad4[] = "QUAD4";
48878569bb4SMatthew G. Knepley       char            elem_type_tri6[] = "TRI6", elem_type_quad9[] = "QUAD9";
48978569bb4SMatthew G. Knepley       char            elem_type_tet4[] = "TET4", elem_type_hex8[] = "HEX8";
49078569bb4SMatthew G. Knepley       char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
49178569bb4SMatthew G. Knepley 
4929566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4939566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4949566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
49578569bb4SMatthew G. Knepley       /* Set Element type */
49678569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
49778569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_tri3;
49878569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tri6;
49978569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
50078569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_quad4;
50178569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_quad9;
50278569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
50378569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_tet4;
50478569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tet10;
50578569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
50678569bb4SMatthew G. Knepley         if (degree == 1) elem_type = elem_type_hex8;
50778569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_hex27;
50878569bb4SMatthew G. Knepley       }
50978569bb4SMatthew G. Knepley       connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
5109566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(PetscMax(27, connectSize) * csSize, &connect));
511792fecdfSBarry Smith       PetscCallExternal(ex_put_block, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1);
51278569bb4SMatthew G. Knepley       /* Find number of vertices, edges, and faces in the closure */
51378569bb4SMatthew G. Knepley       verticesInClosure = nodes[cs][0];
51478569bb4SMatthew G. Knepley       if (depth > 1) {
51578569bb4SMatthew G. Knepley         if (dim == 2) {
5169566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure));
51778569bb4SMatthew G. Knepley         } else if (dim == 3) {
51878569bb4SMatthew G. Knepley           PetscInt *closure = NULL;
51978569bb4SMatthew G. Knepley 
5209566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure));
5219566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
52278569bb4SMatthew G. Knepley           edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
5239566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
52478569bb4SMatthew G. Knepley         }
52578569bb4SMatthew G. Knepley       }
52678569bb4SMatthew G. Knepley       /* Get connectivity for each cell */
52778569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
52878569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
52978569bb4SMatthew G. Knepley         PetscInt  temp, i;
53078569bb4SMatthew G. Knepley 
5319566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
53278569bb4SMatthew G. Knepley         for (i = 0; i < connectSize; ++i) {
53378569bb4SMatthew G. Knepley           if (i < nodes[cs][0]) { /* Vertices */
534fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i + edgesInClosure + facesInClosure + 1) * 2] + 1;
535fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
53678569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1]) { /* Edges */
537fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i - verticesInClosure + facesInClosure + 1) * 2] + 1;
538fe209ef9SBlaise Bourdin             if (nodes[cs][2] == 0) connect[i + off] -= numFaces;
539fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
54078569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3]) { /* Cells */
541fe209ef9SBlaise Bourdin             connect[i + off] = closure[0] + 1;
542fe209ef9SBlaise Bourdin             connect[i + off] -= skipCells;
54378569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0] + nodes[cs][1] + nodes[cs][3] + nodes[cs][2]) { /* Faces */
544fe209ef9SBlaise Bourdin             connect[i + off] = closure[(i - edgesInClosure - verticesInClosure) * 2] + 1;
545fe209ef9SBlaise Bourdin             connect[i + off] -= cellsNotInConnectivity;
54678569bb4SMatthew G. Knepley           } else {
547fe209ef9SBlaise Bourdin             connect[i + off] = -1;
54878569bb4SMatthew G. Knepley           }
54978569bb4SMatthew G. Knepley         }
55078569bb4SMatthew G. Knepley         /* Tetrahedra are inverted */
55178569bb4SMatthew G. Knepley         if (type[cs] == TET) {
5529371c9d4SSatish Balay           temp             = connect[0 + off];
5539371c9d4SSatish Balay           connect[0 + off] = connect[1 + off];
5549371c9d4SSatish Balay           connect[1 + off] = temp;
55578569bb4SMatthew G. Knepley           if (degree == 2) {
5569371c9d4SSatish Balay             temp             = connect[5 + off];
5579371c9d4SSatish Balay             connect[5 + off] = connect[6 + off];
5589371c9d4SSatish Balay             connect[6 + off] = temp;
5599371c9d4SSatish Balay             temp             = connect[7 + off];
5609371c9d4SSatish Balay             connect[7 + off] = connect[8 + off];
5619371c9d4SSatish Balay             connect[8 + off] = temp;
56278569bb4SMatthew G. Knepley           }
56378569bb4SMatthew G. Knepley         }
56478569bb4SMatthew G. Knepley         /* Hexahedra are inverted */
56578569bb4SMatthew G. Knepley         if (type[cs] == HEX) {
5669371c9d4SSatish Balay           temp             = connect[1 + off];
5679371c9d4SSatish Balay           connect[1 + off] = connect[3 + off];
5689371c9d4SSatish Balay           connect[3 + off] = temp;
56978569bb4SMatthew G. Knepley           if (degree == 2) {
5709371c9d4SSatish Balay             temp              = connect[8 + off];
5719371c9d4SSatish Balay             connect[8 + off]  = connect[11 + off];
5729371c9d4SSatish Balay             connect[11 + off] = temp;
5739371c9d4SSatish Balay             temp              = connect[9 + off];
5749371c9d4SSatish Balay             connect[9 + off]  = connect[10 + off];
5759371c9d4SSatish Balay             connect[10 + off] = temp;
5769371c9d4SSatish Balay             temp              = connect[16 + off];
5779371c9d4SSatish Balay             connect[16 + off] = connect[17 + off];
5789371c9d4SSatish Balay             connect[17 + off] = temp;
5799371c9d4SSatish Balay             temp              = connect[18 + off];
5809371c9d4SSatish Balay             connect[18 + off] = connect[19 + off];
5819371c9d4SSatish Balay             connect[19 + off] = temp;
58278569bb4SMatthew G. Knepley 
5839371c9d4SSatish Balay             temp              = connect[12 + off];
5849371c9d4SSatish Balay             connect[12 + off] = connect[16 + off];
5859371c9d4SSatish Balay             connect[16 + off] = temp;
5869371c9d4SSatish Balay             temp              = connect[13 + off];
5879371c9d4SSatish Balay             connect[13 + off] = connect[17 + off];
5889371c9d4SSatish Balay             connect[17 + off] = temp;
5899371c9d4SSatish Balay             temp              = connect[14 + off];
5909371c9d4SSatish Balay             connect[14 + off] = connect[18 + off];
5919371c9d4SSatish Balay             connect[18 + off] = temp;
5929371c9d4SSatish Balay             temp              = connect[15 + off];
5939371c9d4SSatish Balay             connect[15 + off] = connect[19 + off];
5949371c9d4SSatish Balay             connect[19 + off] = temp;
59578569bb4SMatthew G. Knepley 
5969371c9d4SSatish Balay             temp              = connect[23 + off];
5979371c9d4SSatish Balay             connect[23 + off] = connect[26 + off];
5989371c9d4SSatish Balay             connect[26 + off] = temp;
5999371c9d4SSatish Balay             temp              = connect[24 + off];
6009371c9d4SSatish Balay             connect[24 + off] = connect[25 + off];
6019371c9d4SSatish Balay             connect[25 + off] = temp;
6029371c9d4SSatish Balay             temp              = connect[25 + off];
6039371c9d4SSatish Balay             connect[25 + off] = connect[26 + off];
6049371c9d4SSatish Balay             connect[26 + off] = temp;
60578569bb4SMatthew G. Knepley           }
60678569bb4SMatthew G. Knepley         }
607fe209ef9SBlaise Bourdin         off += connectSize;
6089566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
60978569bb4SMatthew G. Knepley       }
610792fecdfSBarry Smith       PetscCallExternal(ex_put_conn, exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0);
61178569bb4SMatthew G. Knepley       skipCells += (nodes[cs][3] == 0) * csSize;
6129566063dSJacob Faibussowitsch       PetscCall(PetscFree(connect));
6139566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
6149566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
61578569bb4SMatthew G. Knepley     }
6169566063dSJacob Faibussowitsch     PetscCall(PetscFree(type));
61778569bb4SMatthew G. Knepley     /* --- Coordinates --- */
6189566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd));
6191e50132fSMatthew G. Knepley     if (num_cs) {
62078569bb4SMatthew G. Knepley       for (d = 0; d < depth; ++d) {
6219566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
62248a46eb9SPierre Jolivet         for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0));
62378569bb4SMatthew G. Knepley       }
6241e50132fSMatthew G. Knepley     }
62578569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
62678569bb4SMatthew G. Knepley       IS              stratumIS;
62778569bb4SMatthew G. Knepley       const PetscInt *cells;
62878569bb4SMatthew G. Knepley       PetscInt        csSize, c;
62978569bb4SMatthew G. Knepley 
6309566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
6319566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
6329566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
63348a46eb9SPierre Jolivet       for (c = 0; c < csSize; ++c) PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0));
6349566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
6359566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
63678569bb4SMatthew G. Knepley     }
637*bdbfaa5dSBlaise Bourdin     if (num_cs) {
6389566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(csIS, &csIdx));
6399566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&csIS));
64078569bb4SMatthew G. Knepley     }
6419566063dSJacob Faibussowitsch     PetscCall(PetscFree(nodes));
6429566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(coordSection));
643*bdbfaa5dSBlaise Bourdin     if (numNodes) {
64478569bb4SMatthew G. Knepley       const char  *coordNames[3] = {"x", "y", "z"};
645233c95e0SFrancesco Ballarin       PetscScalar *closure, *cval;
646233c95e0SFrancesco Ballarin       PetscReal   *coords;
64778569bb4SMatthew G. Knepley       PetscInt     hasDof, n = 0;
64878569bb4SMatthew G. Knepley 
6496823f3c5SBlaise Bourdin       /* There can't be more than 24 values in the closure of a point for the coord coordSection */
6509566063dSJacob Faibussowitsch       PetscCall(PetscCalloc3(numNodes * 3, &coords, dim, &cval, 24, &closure));
6519566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord));
6529566063dSJacob Faibussowitsch       PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
65378569bb4SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
6549566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(coordSection, p, &hasDof));
65578569bb4SMatthew G. Knepley         if (hasDof) {
65678569bb4SMatthew G. Knepley           PetscInt closureSize = 24, j;
65778569bb4SMatthew G. Knepley 
6589566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure));
65978569bb4SMatthew G. Knepley           for (d = 0; d < dim; ++d) {
66078569bb4SMatthew G. Knepley             cval[d] = 0.0;
66178569bb4SMatthew G. Knepley             for (j = 0; j < closureSize / dim; j++) cval[d] += closure[j * dim + d];
662233c95e0SFrancesco Ballarin             coords[d * numNodes + n] = PetscRealPart(cval[d]) * dim / closureSize;
66378569bb4SMatthew G. Knepley           }
66478569bb4SMatthew G. Knepley           ++n;
66578569bb4SMatthew G. Knepley         }
66678569bb4SMatthew G. Knepley       }
667792fecdfSBarry Smith       PetscCallExternal(ex_put_coord, exo->exoid, &coords[0 * numNodes], &coords[1 * numNodes], &coords[2 * numNodes]);
6689566063dSJacob Faibussowitsch       PetscCall(PetscFree3(coords, cval, closure));
669792fecdfSBarry Smith       PetscCallExternal(ex_put_coord_names, exo->exoid, (char **)coordNames);
67078569bb4SMatthew G. Knepley     }
6716823f3c5SBlaise Bourdin 
672fe209ef9SBlaise Bourdin     /* --- Node Sets/Vertex Sets --- */
673fe209ef9SBlaise Bourdin     DMHasLabel(dm, "Vertex Sets", &hasLabel);
674fe209ef9SBlaise Bourdin     if (hasLabel) {
675fe209ef9SBlaise Bourdin       PetscInt        i, vs, vsSize;
676fe209ef9SBlaise Bourdin       const PetscInt *vsIdx, *vertices;
677fe209ef9SBlaise Bourdin       PetscInt       *nodeList;
678fe209ef9SBlaise Bourdin       IS              vsIS, stratumIS;
679fe209ef9SBlaise Bourdin       DMLabel         vsLabel;
6809566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel));
6819566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(vsLabel, &vsIS));
6829566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(vsIS, &vsIdx));
683fe209ef9SBlaise Bourdin       for (vs = 0; vs < num_vs; ++vs) {
6849566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS));
6859566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(stratumIS, &vertices));
6869566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &vsSize));
6879566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(vsSize, &nodeList));
688ad540459SPierre Jolivet         for (i = 0; i < vsSize; ++i) nodeList[i] = vertices[i] - skipCells + 1;
689792fecdfSBarry Smith         PetscCallExternal(ex_put_set_param, exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0);
690792fecdfSBarry Smith         PetscCallExternal(ex_put_set, exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL);
6919566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(stratumIS, &vertices));
6929566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
6939566063dSJacob Faibussowitsch         PetscCall(PetscFree(nodeList));
694fe209ef9SBlaise Bourdin       }
6959566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(vsIS, &vsIdx));
6969566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&vsIS));
697fe209ef9SBlaise Bourdin     }
698fe209ef9SBlaise Bourdin     /* --- Side Sets/Face Sets --- */
6999566063dSJacob Faibussowitsch     PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel));
700fe209ef9SBlaise Bourdin     if (hasLabel) {
701fe209ef9SBlaise Bourdin       PetscInt        i, j, fs, fsSize;
702fe209ef9SBlaise Bourdin       const PetscInt *fsIdx, *faces;
703fe209ef9SBlaise Bourdin       IS              fsIS, stratumIS;
704fe209ef9SBlaise Bourdin       DMLabel         fsLabel;
705fe209ef9SBlaise Bourdin       PetscInt        numPoints, *points;
706fe209ef9SBlaise Bourdin       PetscInt        elem_list_size = 0;
707fe209ef9SBlaise Bourdin       PetscInt       *elem_list, *elem_ind, *side_list;
708fe209ef9SBlaise Bourdin 
7099566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel));
710fe209ef9SBlaise Bourdin       /* Compute size of Node List and Element List */
7119566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(fsLabel, &fsIS));
7129566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fsIS, &fsIdx));
713fe209ef9SBlaise Bourdin       for (fs = 0; fs < num_fs; ++fs) {
7149566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
7159566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &fsSize));
716fe209ef9SBlaise Bourdin         elem_list_size += fsSize;
7179566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
718fe209ef9SBlaise Bourdin       }
719*bdbfaa5dSBlaise Bourdin       if (num_fs) {
720d0609cedSBarry Smith         PetscCall(PetscMalloc3(num_fs, &elem_ind, elem_list_size, &elem_list, elem_list_size, &side_list));
721fe209ef9SBlaise Bourdin         elem_ind[0] = 0;
722fe209ef9SBlaise Bourdin         for (fs = 0; fs < num_fs; ++fs) {
7239566063dSJacob Faibussowitsch           PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
7249566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(stratumIS, &faces));
7259566063dSJacob Faibussowitsch           PetscCall(ISGetSize(stratumIS, &fsSize));
726fe209ef9SBlaise Bourdin           /* Set Parameters */
727792fecdfSBarry Smith           PetscCallExternal(ex_put_set_param, exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0);
728fe209ef9SBlaise Bourdin           /* Indices */
729ad540459SPierre Jolivet           if (fs < num_fs - 1) elem_ind[fs + 1] = elem_ind[fs] + fsSize;
730fe209ef9SBlaise Bourdin 
731fe209ef9SBlaise Bourdin           for (i = 0; i < fsSize; ++i) {
732fe209ef9SBlaise Bourdin             /* Element List */
733fe209ef9SBlaise Bourdin             points = NULL;
7349566063dSJacob Faibussowitsch             PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
735fe209ef9SBlaise Bourdin             elem_list[elem_ind[fs] + i] = points[2] + 1;
7369566063dSJacob Faibussowitsch             PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
737fe209ef9SBlaise Bourdin 
738fe209ef9SBlaise Bourdin             /* Side List */
739fe209ef9SBlaise Bourdin             points = NULL;
7409566063dSJacob Faibussowitsch             PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
741fe209ef9SBlaise Bourdin             for (j = 1; j < numPoints; ++j) {
742ad540459SPierre Jolivet               if (points[j * 2] == faces[i]) break;
743fe209ef9SBlaise Bourdin             }
744fe209ef9SBlaise Bourdin             /* Convert HEX sides */
745fe209ef9SBlaise Bourdin             if (numPoints == 27) {
7469371c9d4SSatish Balay               if (j == 1) {
7479371c9d4SSatish Balay                 j = 5;
7489371c9d4SSatish Balay               } else if (j == 2) {
7499371c9d4SSatish Balay                 j = 6;
7509371c9d4SSatish Balay               } else if (j == 3) {
7519371c9d4SSatish Balay                 j = 1;
7529371c9d4SSatish Balay               } else if (j == 4) {
7539371c9d4SSatish Balay                 j = 3;
7549371c9d4SSatish Balay               } else if (j == 5) {
7559371c9d4SSatish Balay                 j = 2;
7569371c9d4SSatish Balay               } else if (j == 6) {
7579371c9d4SSatish Balay                 j = 4;
7589371c9d4SSatish Balay               }
759fe209ef9SBlaise Bourdin             }
760fe209ef9SBlaise Bourdin             /* Convert TET sides */
761fe209ef9SBlaise Bourdin             if (numPoints == 15) {
762fe209ef9SBlaise Bourdin               --j;
763ad540459SPierre Jolivet               if (j == 0) j = 4;
764fe209ef9SBlaise Bourdin             }
765fe209ef9SBlaise Bourdin             side_list[elem_ind[fs] + i] = j;
7669566063dSJacob Faibussowitsch             PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i] - 1, PETSC_TRUE, &numPoints, &points));
767fe209ef9SBlaise Bourdin           }
7689566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(stratumIS, &faces));
7699566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&stratumIS));
770fe209ef9SBlaise Bourdin         }
7719566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(fsIS, &fsIdx));
7729566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&fsIS));
773fe209ef9SBlaise Bourdin 
774fe209ef9SBlaise Bourdin         /* Put side sets */
77548a46eb9SPierre 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]]);
7769566063dSJacob Faibussowitsch         PetscCall(PetscFree3(elem_ind, elem_list, side_list));
777fe209ef9SBlaise Bourdin       }
778*bdbfaa5dSBlaise Bourdin     }
7796823f3c5SBlaise Bourdin     /*
7806823f3c5SBlaise Bourdin       close the exodus file
7816823f3c5SBlaise Bourdin     */
7826823f3c5SBlaise Bourdin     ex_close(exo->exoid);
7836823f3c5SBlaise Bourdin     exo->exoid = -1;
7846823f3c5SBlaise Bourdin   }
7859566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&coordSection));
7866823f3c5SBlaise Bourdin 
7876823f3c5SBlaise Bourdin   /*
7886823f3c5SBlaise Bourdin     reopen the file in parallel
7896823f3c5SBlaise Bourdin   */
7906823f3c5SBlaise Bourdin   EXO_mode = EX_WRITE;
7916823f3c5SBlaise Bourdin   #if defined(PETSC_USE_64BIT_INDICES)
7926823f3c5SBlaise Bourdin   EXO_mode += EX_ALL_INT64_API;
7936823f3c5SBlaise Bourdin   #endif
7946823f3c5SBlaise Bourdin   CPU_word_size = sizeof(PetscReal);
7956823f3c5SBlaise Bourdin   IO_word_size  = sizeof(PetscReal);
7966823f3c5SBlaise Bourdin   exo->exoid    = ex_open_par(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size, &EXO_version, comm, MPI_INFO_NULL);
79708401ef6SPierre Jolivet   PetscCheck(exo->exoid >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
7986823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
7996823f3c5SBlaise Bourdin }
8006823f3c5SBlaise Bourdin 
8016823f3c5SBlaise Bourdin /*
8026823f3c5SBlaise Bourdin   VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
8036823f3c5SBlaise Bourdin 
8046823f3c5SBlaise Bourdin   Collective on v
8056823f3c5SBlaise Bourdin 
8066823f3c5SBlaise Bourdin   Input Parameters:
8076823f3c5SBlaise Bourdin + v  - The vector to be written
8086823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
8096823f3c5SBlaise Bourdin 
8106823f3c5SBlaise Bourdin   Notes:
8116823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
8126823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
8136823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
8146823f3c5SBlaise Bourdin   amongst all variables.
8156823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
8166823f3c5SBlaise Bourdin   possibly corrupting the file
8176823f3c5SBlaise Bourdin 
8186823f3c5SBlaise Bourdin   Level: beginner
8196823f3c5SBlaise Bourdin 
820c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
8216823f3c5SBlaise Bourdin @*/
822d71ae5a4SJacob Faibussowitsch PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer)
823d71ae5a4SJacob Faibussowitsch {
8246823f3c5SBlaise Bourdin   DM          dm;
8256823f3c5SBlaise Bourdin   MPI_Comm    comm;
8266823f3c5SBlaise Bourdin   PetscMPIInt rank;
8276823f3c5SBlaise Bourdin 
8286823f3c5SBlaise Bourdin   int         exoid, offsetN = 0, offsetZ = 0;
8296823f3c5SBlaise Bourdin   const char *vecname;
8306823f3c5SBlaise Bourdin   PetscInt    step;
8316823f3c5SBlaise Bourdin 
8326823f3c5SBlaise Bourdin   PetscFunctionBegin;
8339566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
8349566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8359566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
8369566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8379566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
8386823f3c5SBlaise Bourdin 
8399566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
8409566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN));
8419566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ));
8421dca8a05SBarry Smith   PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
8436823f3c5SBlaise Bourdin   if (offsetN > 0) {
8449566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN));
8456823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
8469566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ));
84798921bdaSJacob Faibussowitsch   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
8486823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
8496823f3c5SBlaise Bourdin }
8506823f3c5SBlaise Bourdin 
8516823f3c5SBlaise Bourdin /*
8526823f3c5SBlaise Bourdin   VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
8536823f3c5SBlaise Bourdin 
8546823f3c5SBlaise Bourdin   Collective on v
8556823f3c5SBlaise Bourdin 
8566823f3c5SBlaise Bourdin   Input Parameters:
8576823f3c5SBlaise Bourdin + v  - The vector to be written
8586823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
8596823f3c5SBlaise Bourdin 
8606823f3c5SBlaise Bourdin   Notes:
8616823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
8626823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
8636823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
8646823f3c5SBlaise Bourdin   amongst all variables.
8656823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
8666823f3c5SBlaise Bourdin   possibly corrupting the file
8676823f3c5SBlaise Bourdin 
8686823f3c5SBlaise Bourdin   Level: beginner
8696823f3c5SBlaise Bourdin 
870c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII()`, `VecView_PlexExodusII()`
8716823f3c5SBlaise Bourdin @*/
872d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer)
873d71ae5a4SJacob Faibussowitsch {
8746823f3c5SBlaise Bourdin   DM          dm;
8756823f3c5SBlaise Bourdin   MPI_Comm    comm;
8766823f3c5SBlaise Bourdin   PetscMPIInt rank;
8776823f3c5SBlaise Bourdin 
8786823f3c5SBlaise Bourdin   int         exoid, offsetN = 0, offsetZ = 0;
8796823f3c5SBlaise Bourdin   const char *vecname;
8806823f3c5SBlaise Bourdin   PetscInt    step;
8816823f3c5SBlaise Bourdin 
8826823f3c5SBlaise Bourdin   PetscFunctionBegin;
8839566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
8849566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8859566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer, &exoid));
8869566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8879566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject)v, &vecname));
8886823f3c5SBlaise Bourdin 
8899566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm, &step, NULL));
8909566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_NODAL, vecname, &offsetN));
8919566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid, EX_ELEM_BLOCK, vecname, &offsetZ));
8921dca8a05SBarry Smith   PetscCheck(offsetN > 0 || offsetZ > 0, comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
8931dca8a05SBarry Smith   if (offsetN > 0) PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v, exoid, (int)step + 1, offsetN));
8941dca8a05SBarry Smith   else if (offsetZ > 0) PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v, exoid, (int)step + 1, offsetZ));
8951dca8a05SBarry Smith   else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
89678569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
89778569bb4SMatthew G. Knepley }
89878569bb4SMatthew G. Knepley 
89978569bb4SMatthew G. Knepley /*
90078569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
90178569bb4SMatthew G. Knepley 
90278569bb4SMatthew G. Knepley   Collective on v
90378569bb4SMatthew G. Knepley 
90478569bb4SMatthew G. Knepley   Input Parameters:
90578569bb4SMatthew G. Knepley + v  - The vector to be written
90678569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
9076823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
9086823f3c5SBlaise Bourdin - offset - the location of the variable in the file
90978569bb4SMatthew G. Knepley 
91078569bb4SMatthew G. Knepley   Notes:
91178569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
91278569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
91378569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
91478569bb4SMatthew G. Knepley   amongst all nodal variables.
91578569bb4SMatthew G. Knepley 
91678569bb4SMatthew G. Knepley   Level: beginner
91778569bb4SMatthew G. Knepley 
918c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecLoadNodal_PlexEXO()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
91978569bb4SMatthew G. Knepley @*/
920d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
921d71ae5a4SJacob Faibussowitsch {
92278569bb4SMatthew G. Knepley   MPI_Comm           comm;
92378569bb4SMatthew G. Knepley   PetscMPIInt        size;
92478569bb4SMatthew G. Knepley   DM                 dm;
92578569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
92622a7544eSVaclav Hapla   const PetscScalar *varray;
92778569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
92878569bb4SMatthew G. Knepley   PetscBool          useNatural;
92978569bb4SMatthew G. Knepley 
93078569bb4SMatthew G. Knepley   PetscFunctionBegin;
9319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
9329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9339566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
9349566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
93578569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
93678569bb4SMatthew G. Knepley   if (useNatural) {
937b7352c5cSAlexis Marboeuf     PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
9389566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
9399566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
94078569bb4SMatthew G. Knepley   } else {
94178569bb4SMatthew G. Knepley     vNatural = v;
94278569bb4SMatthew G. Knepley   }
9436823f3c5SBlaise Bourdin 
94478569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
94578569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
94678569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
9479566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
9489566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
94978569bb4SMatthew G. Knepley   if (bs == 1) {
9509566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(vNatural, &varray));
951792fecdfSBarry Smith     PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
9529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(vNatural, &varray));
95378569bb4SMatthew G. Knepley   } else {
95478569bb4SMatthew G. Knepley     IS       compIS;
95578569bb4SMatthew G. Knepley     PetscInt c;
95678569bb4SMatthew G. Knepley 
9579566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
95878569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
9599566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
9609566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
9619566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vComp, &varray));
962792fecdfSBarry Smith       PetscCallExternal(ex_put_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
9639566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vComp, &varray));
9649566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
96578569bb4SMatthew G. Knepley     }
9669566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
96778569bb4SMatthew G. Knepley   }
968b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(VecDestroy(&vNatural));
96978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
97078569bb4SMatthew G. Knepley }
97178569bb4SMatthew G. Knepley 
97278569bb4SMatthew G. Knepley /*
97378569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
97478569bb4SMatthew G. Knepley 
97578569bb4SMatthew G. Knepley   Collective on v
97678569bb4SMatthew G. Knepley 
97778569bb4SMatthew G. Knepley   Input Parameters:
97878569bb4SMatthew G. Knepley + v  - The vector to be written
97978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
9806823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
9816823f3c5SBlaise Bourdin - offset - the location of the variable in the file
98278569bb4SMatthew G. Knepley 
98378569bb4SMatthew G. Knepley   Notes:
98478569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
98578569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
98678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
98778569bb4SMatthew G. Knepley   amongst all nodal variables.
98878569bb4SMatthew G. Knepley 
98978569bb4SMatthew G. Knepley   Level: beginner
99078569bb4SMatthew G. Knepley 
991db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecViewZonal_PlexEXO()`, `VecLoadZonal_PlexEXO()`
99278569bb4SMatthew G. Knepley */
993d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
994d71ae5a4SJacob Faibussowitsch {
99578569bb4SMatthew G. Knepley   MPI_Comm     comm;
99678569bb4SMatthew G. Knepley   PetscMPIInt  size;
99778569bb4SMatthew G. Knepley   DM           dm;
99878569bb4SMatthew G. Knepley   Vec          vNatural, vComp;
99922a7544eSVaclav Hapla   PetscScalar *varray;
100078569bb4SMatthew G. Knepley   PetscInt     xs, xe, bs;
100178569bb4SMatthew G. Knepley   PetscBool    useNatural;
100278569bb4SMatthew G. Knepley 
100378569bb4SMatthew G. Knepley   PetscFunctionBegin;
10049566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
10059566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10069566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
10079566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
100878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
1009b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
1010ad540459SPierre Jolivet   else vNatural = v;
10116823f3c5SBlaise Bourdin 
101278569bb4SMatthew G. Knepley   /* Read local chunk from the file */
10139566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
10149566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
101578569bb4SMatthew G. Knepley   if (bs == 1) {
10169566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vNatural, &varray));
1017792fecdfSBarry Smith     PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset, 1, xs + 1, xe - xs, varray);
10189566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vNatural, &varray));
101978569bb4SMatthew G. Knepley   } else {
102078569bb4SMatthew G. Knepley     IS       compIS;
102178569bb4SMatthew G. Knepley     PetscInt c;
102278569bb4SMatthew G. Knepley 
10239566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
102478569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
10259566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
10269566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
10279566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vComp, &varray));
1028792fecdfSBarry Smith       PetscCallExternal(ex_get_partial_var, exoid, step, EX_NODAL, offset + c, 1, xs / bs + 1, (xe - xs) / bs, varray);
10299566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vComp, &varray));
10309566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
103178569bb4SMatthew G. Knepley     }
10329566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
103378569bb4SMatthew G. Knepley   }
103478569bb4SMatthew G. Knepley   if (useNatural) {
10359566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
10369566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
1037b7352c5cSAlexis Marboeuf     PetscCall(VecDestroy(&vNatural));
103878569bb4SMatthew G. Knepley   }
103978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
104078569bb4SMatthew G. Knepley }
104178569bb4SMatthew G. Knepley 
104278569bb4SMatthew G. Knepley /*
104378569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
104478569bb4SMatthew G. Knepley 
104578569bb4SMatthew G. Knepley   Collective on v
104678569bb4SMatthew G. Knepley 
104778569bb4SMatthew G. Knepley   Input Parameters:
104878569bb4SMatthew G. Knepley + v  - The vector to be written
104978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
10506823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
10516823f3c5SBlaise Bourdin - offset - the location of the variable in the file
105278569bb4SMatthew G. Knepley 
105378569bb4SMatthew G. Knepley   Notes:
105478569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
105578569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
105678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
105778569bb4SMatthew G. Knepley   amongst all zonal variables.
105878569bb4SMatthew G. Knepley 
105978569bb4SMatthew G. Knepley   Level: beginner
106078569bb4SMatthew G. Knepley 
1061c2e3fba1SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
106278569bb4SMatthew G. Knepley */
1063d71ae5a4SJacob Faibussowitsch PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
1064d71ae5a4SJacob Faibussowitsch {
106578569bb4SMatthew G. Knepley   MPI_Comm           comm;
106678569bb4SMatthew G. Knepley   PetscMPIInt        size;
106778569bb4SMatthew G. Knepley   DM                 dm;
106878569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
106922a7544eSVaclav Hapla   const PetscScalar *varray;
107078569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
107178569bb4SMatthew G. Knepley   PetscBool          useNatural;
107278569bb4SMatthew G. Knepley   IS                 compIS;
107378569bb4SMatthew G. Knepley   PetscInt          *csSize, *csID;
107478569bb4SMatthew G. Knepley   PetscInt           numCS, set, csxs = 0;
107578569bb4SMatthew G. Knepley 
107678569bb4SMatthew G. Knepley   PetscFunctionBegin;
10779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
10789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10799566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
10809566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
108178569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
108278569bb4SMatthew G. Knepley   if (useNatural) {
1083b7352c5cSAlexis Marboeuf     PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
10849566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
10859566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
108678569bb4SMatthew G. Knepley   } else {
108778569bb4SMatthew G. Knepley     vNatural = v;
108878569bb4SMatthew G. Knepley   }
10896823f3c5SBlaise Bourdin 
109078569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
109178569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
109278569bb4SMatthew G. Knepley      We assume that they are stored sequentially
109378569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1094a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
109578569bb4SMatthew G. Knepley      to figure out what to save where. */
109678569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
10979566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1098792fecdfSBarry Smith   PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
109978569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
110078569bb4SMatthew G. Knepley     ex_block block;
110178569bb4SMatthew G. Knepley 
110278569bb4SMatthew G. Knepley     block.id   = csID[set];
110378569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1104792fecdfSBarry Smith     PetscCallExternal(ex_get_block_param, exoid, &block);
110578569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
110678569bb4SMatthew G. Knepley   }
11079566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
11089566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
11099566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
111078569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
111178569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
111278569bb4SMatthew G. Knepley 
111378569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
111478569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
111578569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
111678569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
111778569bb4SMatthew G. Knepley     if (bs == 1) {
11189566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vNatural, &varray));
1119792fecdfSBarry 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)]);
11209566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vNatural, &varray));
112178569bb4SMatthew G. Knepley     } else {
112278569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
11239566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
11249566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
11259566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vComp, &varray));
1126792fecdfSBarry 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)]);
11279566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vComp, &varray));
11289566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
112978569bb4SMatthew G. Knepley       }
113078569bb4SMatthew G. Knepley     }
113178569bb4SMatthew G. Knepley     csxs += csSize[set];
113278569bb4SMatthew G. Knepley   }
11339566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
11349566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
1135b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(VecDestroy(&vNatural));
113678569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
113778569bb4SMatthew G. Knepley }
113878569bb4SMatthew G. Knepley 
113978569bb4SMatthew G. Knepley /*
114078569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
114178569bb4SMatthew G. Knepley 
114278569bb4SMatthew G. Knepley   Collective on v
114378569bb4SMatthew G. Knepley 
114478569bb4SMatthew G. Knepley   Input Parameters:
114578569bb4SMatthew G. Knepley + v  - The vector to be written
114678569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
11476823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
11486823f3c5SBlaise Bourdin - offset - the location of the variable in the file
114978569bb4SMatthew G. Knepley 
115078569bb4SMatthew G. Knepley   Notes:
115178569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
115278569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
115378569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
115478569bb4SMatthew G. Knepley   amongst all zonal variables.
115578569bb4SMatthew G. Knepley 
115678569bb4SMatthew G. Knepley   Level: beginner
115778569bb4SMatthew G. Knepley 
1158db781477SPatrick Sanan .seealso: `EXOGetVarIndex_Internal()`, `DMPlexView_ExodusII_Internal()`, `VecViewPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Nodal_Internal()`, `VecLoadPlex_ExodusII_Zonal_Internal()`
115978569bb4SMatthew G. Knepley */
1160d71ae5a4SJacob Faibussowitsch PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
1161d71ae5a4SJacob Faibussowitsch {
116278569bb4SMatthew G. Knepley   MPI_Comm     comm;
116378569bb4SMatthew G. Knepley   PetscMPIInt  size;
116478569bb4SMatthew G. Knepley   DM           dm;
116578569bb4SMatthew G. Knepley   Vec          vNatural, vComp;
116622a7544eSVaclav Hapla   PetscScalar *varray;
116778569bb4SMatthew G. Knepley   PetscInt     xs, xe, bs;
116878569bb4SMatthew G. Knepley   PetscBool    useNatural;
116978569bb4SMatthew G. Knepley   IS           compIS;
117078569bb4SMatthew G. Knepley   PetscInt    *csSize, *csID;
117178569bb4SMatthew G. Knepley   PetscInt     numCS, set, csxs = 0;
117278569bb4SMatthew G. Knepley 
117378569bb4SMatthew G. Knepley   PetscFunctionBegin;
11749566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
11759566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
11769566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
11779566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
117878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
1179b7352c5cSAlexis Marboeuf   if (useNatural) PetscCall(DMPlexCreateNaturalVector(dm, &vNatural));
1180ad540459SPierre Jolivet   else vNatural = v;
11816823f3c5SBlaise Bourdin 
118278569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
118378569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
118478569bb4SMatthew G. Knepley      We assume that they are stored sequentially
118578569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1186a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
118778569bb4SMatthew G. Knepley      to figure out what to save where. */
118878569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
11899566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1190792fecdfSBarry Smith   PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, csID);
119178569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
119278569bb4SMatthew G. Knepley     ex_block block;
119378569bb4SMatthew G. Knepley 
119478569bb4SMatthew G. Knepley     block.id   = csID[set];
119578569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1196792fecdfSBarry Smith     PetscCallExternal(ex_get_block_param, exoid, &block);
119778569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
119878569bb4SMatthew G. Knepley   }
11999566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
12009566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
12019566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe - xs) / bs, xs, bs, &compIS));
120278569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
120378569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
120478569bb4SMatthew G. Knepley 
120578569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
120678569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
120778569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
120878569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe / bs, csxs + csSize[set]) - PetscMax(xs / bs, csxs));
120978569bb4SMatthew G. Knepley     if (bs == 1) {
12109566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vNatural, &varray));
1211792fecdfSBarry 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)]);
12129566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vNatural, &varray));
121378569bb4SMatthew G. Knepley     } else {
121478569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
12159566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe - xs) / bs, xs + c, bs));
12169566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
12179566063dSJacob Faibussowitsch         PetscCall(VecGetArray(vComp, &varray));
1218792fecdfSBarry 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)]);
12199566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(vComp, &varray));
12209566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
122178569bb4SMatthew G. Knepley       }
122278569bb4SMatthew G. Knepley     }
122378569bb4SMatthew G. Knepley     csxs += csSize[set];
122478569bb4SMatthew G. Knepley   }
12259566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
12269566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
122778569bb4SMatthew G. Knepley   if (useNatural) {
12289566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
12299566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
1230b7352c5cSAlexis Marboeuf     PetscCall(VecDestroy(&vNatural));
123178569bb4SMatthew G. Knepley   }
123278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
123378569bb4SMatthew G. Knepley }
1234b53c8456SSatish Balay #endif
123578569bb4SMatthew G. Knepley 
12361e50132fSMatthew G. Knepley /*@
1237a1cb98faSBarry Smith   PetscViewerExodusIIGetId - Get the file id of the `PETSCVIEWEREXODUSII` file
12381e50132fSMatthew G. Knepley 
1239a1cb98faSBarry Smith   Logically Collective on viewer
12401e50132fSMatthew G. Knepley 
12411e50132fSMatthew G. Knepley   Input Parameter:
1242a1cb98faSBarry Smith .  viewer - the `PetscViewer`
12431e50132fSMatthew G. Knepley 
12441e50132fSMatthew G. Knepley   Output Parameter:
1245d8d19677SJose E. Roman .  exoid - The ExodusII file id
12461e50132fSMatthew G. Knepley 
12471e50132fSMatthew G. Knepley   Level: intermediate
12481e50132fSMatthew G. Knepley 
1249a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
12501e50132fSMatthew G. Knepley @*/
1251d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid)
1252d71ae5a4SJacob Faibussowitsch {
12531e50132fSMatthew G. Knepley   PetscFunctionBegin;
12541e50132fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1255cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetId_C", (PetscViewer, int *), (viewer, exoid));
12566823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12576823f3c5SBlaise Bourdin }
12586823f3c5SBlaise Bourdin 
12596823f3c5SBlaise Bourdin /*@
12606823f3c5SBlaise Bourdin    PetscViewerExodusIISetOrder - Set the elements order in the exodusII file.
12616823f3c5SBlaise Bourdin 
12626823f3c5SBlaise Bourdin    Collective
12636823f3c5SBlaise Bourdin 
12646823f3c5SBlaise Bourdin    Input Parameters:
1265a1cb98faSBarry Smith +  viewer - the `PETSCVIEWEREXODUSII` viewer
126698a6a78aSPierre Jolivet -  order - elements order
12676823f3c5SBlaise Bourdin 
12686823f3c5SBlaise Bourdin    Output Parameter:
12696823f3c5SBlaise Bourdin 
12706823f3c5SBlaise Bourdin    Level: beginner
12716823f3c5SBlaise Bourdin 
1272a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
12736823f3c5SBlaise Bourdin @*/
1274d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order)
1275d71ae5a4SJacob Faibussowitsch {
12766823f3c5SBlaise Bourdin   PetscFunctionBegin;
12776823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1278cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerSetOrder_C", (PetscViewer, PetscInt), (viewer, order));
12796823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12806823f3c5SBlaise Bourdin }
12816823f3c5SBlaise Bourdin 
12826823f3c5SBlaise Bourdin /*@
12836823f3c5SBlaise Bourdin    PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file.
12846823f3c5SBlaise Bourdin 
12856823f3c5SBlaise Bourdin    Collective
12866823f3c5SBlaise Bourdin 
12876823f3c5SBlaise Bourdin    Input Parameters:
1288a1cb98faSBarry Smith +  viewer - the `PETSCVIEWEREXODUSII` viewer
128998a6a78aSPierre Jolivet -  order - elements order
12906823f3c5SBlaise Bourdin 
12916823f3c5SBlaise Bourdin    Output Parameter:
12926823f3c5SBlaise Bourdin 
12936823f3c5SBlaise Bourdin    Level: beginner
12946823f3c5SBlaise Bourdin 
1295a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerExodusIIGetId()`, `PetscViewerExodusIIGetOrder()`, `PetscViewerExodusIISetOrder()`
12966823f3c5SBlaise Bourdin @*/
1297d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order)
1298d71ae5a4SJacob Faibussowitsch {
12996823f3c5SBlaise Bourdin   PetscFunctionBegin;
13006823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1301cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetOrder_C", (PetscViewer, PetscInt *), (viewer, order));
13026823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
13036823f3c5SBlaise Bourdin }
13046823f3c5SBlaise Bourdin 
13056823f3c5SBlaise Bourdin /*@C
13066823f3c5SBlaise Bourdin    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
13076823f3c5SBlaise Bourdin 
13086823f3c5SBlaise Bourdin    Collective
13096823f3c5SBlaise Bourdin 
13106823f3c5SBlaise Bourdin    Input Parameters:
13116823f3c5SBlaise Bourdin +  comm - MPI communicator
13126823f3c5SBlaise Bourdin .  name - name of file
13136823f3c5SBlaise Bourdin -  type - type of file
1314a1cb98faSBarry Smith .vb
1315a1cb98faSBarry Smith     FILE_MODE_WRITE - create new file for binary output
1316a1cb98faSBarry Smith     FILE_MODE_READ - open existing file for binary input
1317a1cb98faSBarry Smith     FILE_MODE_APPEND - open existing file for binary output
1318a1cb98faSBarry Smith .ve
13196823f3c5SBlaise Bourdin 
13206823f3c5SBlaise Bourdin    Output Parameter:
1321a1cb98faSBarry Smith .  exo - `PETSCVIEWEREXODUSII` `PetscViewer` for Exodus II input/output to use with the specified file
13226823f3c5SBlaise Bourdin 
13236823f3c5SBlaise Bourdin    Level: beginner
13246823f3c5SBlaise Bourdin 
1325a1cb98faSBarry Smith .seealso: `PETSCVIEWEREXODUSII`, `PetscViewer`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
1326db781477SPatrick Sanan           `DMLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
13276823f3c5SBlaise Bourdin @*/
1328d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo)
1329d71ae5a4SJacob Faibussowitsch {
13306823f3c5SBlaise Bourdin   PetscFunctionBegin;
13319566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, exo));
13329566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII));
13339566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*exo, type));
13349566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*exo, name));
13359566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*exo));
13361e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
13371e50132fSMatthew G. Knepley }
13381e50132fSMatthew G. Knepley 
133933751fbdSMichael Lange /*@C
1340a1cb98faSBarry Smith   DMPlexCreateExodusFromFile - Create a `DMPLEX` mesh from an ExodusII file.
134133751fbdSMichael Lange 
1342d083f849SBarry Smith   Collective
134333751fbdSMichael Lange 
134433751fbdSMichael Lange   Input Parameters:
134533751fbdSMichael Lange + comm  - The MPI communicator
134633751fbdSMichael Lange . filename - The name of the ExodusII file
134733751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
134833751fbdSMichael Lange 
134933751fbdSMichael Lange   Output Parameter:
1350a1cb98faSBarry Smith . dm  - The `DM` object representing the mesh
135133751fbdSMichael Lange 
135233751fbdSMichael Lange   Level: beginner
135333751fbdSMichael Lange 
1354a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMCreate()`, `DMPlexCreateExodus()`
135533751fbdSMichael Lange @*/
1356d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
1357d71ae5a4SJacob Faibussowitsch {
135833751fbdSMichael Lange   PetscMPIInt rank;
135933751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1360e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
136133751fbdSMichael Lange   float version;
136233751fbdSMichael Lange #endif
136333751fbdSMichael Lange 
136433751fbdSMichael Lange   PetscFunctionBegin;
136533751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
13669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
136733751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1368dd400576SPatrick Sanan   if (rank == 0) {
136933751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
137008401ef6SPierre Jolivet     PetscCheck(exoid > 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
137133751fbdSMichael Lange   }
13729566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm));
137348a46eb9SPierre Jolivet   if (rank == 0) PetscCallExternal(ex_close, exoid);
1374b458e8f1SJose E. Roman   PetscFunctionReturn(0);
137533751fbdSMichael Lange #else
137633751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
137733751fbdSMichael Lange #endif
137833751fbdSMichael Lange }
137933751fbdSMichael Lange 
13808f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
1381d71ae5a4SJacob Faibussowitsch static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct)
1382d71ae5a4SJacob Faibussowitsch {
13838f861fbcSMatthew G. Knepley   PetscBool flg;
13848f861fbcSMatthew G. Knepley 
13858f861fbcSMatthew G. Knepley   PetscFunctionBegin;
13868f861fbcSMatthew G. Knepley   *ct = DM_POLYTOPE_UNKNOWN;
13879566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI", &flg));
13889371c9d4SSatish Balay   if (flg) {
13899371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRIANGLE;
13909371c9d4SSatish Balay     goto done;
13919371c9d4SSatish Balay   }
13929566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI3", &flg));
13939371c9d4SSatish Balay   if (flg) {
13949371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRIANGLE;
13959371c9d4SSatish Balay     goto done;
13969371c9d4SSatish Balay   }
13979566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD", &flg));
13989371c9d4SSatish Balay   if (flg) {
13999371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
14009371c9d4SSatish Balay     goto done;
14019371c9d4SSatish Balay   }
14029566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg));
14039371c9d4SSatish Balay   if (flg) {
14049371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
14059371c9d4SSatish Balay     goto done;
14069371c9d4SSatish Balay   }
14079566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg));
14089371c9d4SSatish Balay   if (flg) {
14099371c9d4SSatish Balay     *ct = DM_POLYTOPE_QUADRILATERAL;
14109371c9d4SSatish Balay     goto done;
14119371c9d4SSatish Balay   }
14129566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TETRA", &flg));
14139371c9d4SSatish Balay   if (flg) {
14149371c9d4SSatish Balay     *ct = DM_POLYTOPE_TETRAHEDRON;
14159371c9d4SSatish Balay     goto done;
14169371c9d4SSatish Balay   }
14179566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TET4", &flg));
14189371c9d4SSatish Balay   if (flg) {
14199371c9d4SSatish Balay     *ct = DM_POLYTOPE_TETRAHEDRON;
14209371c9d4SSatish Balay     goto done;
14219371c9d4SSatish Balay   }
14229566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg));
14239371c9d4SSatish Balay   if (flg) {
14249371c9d4SSatish Balay     *ct = DM_POLYTOPE_TRI_PRISM;
14259371c9d4SSatish Balay     goto done;
14269371c9d4SSatish Balay   }
14279566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX", &flg));
14289371c9d4SSatish Balay   if (flg) {
14299371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
14309371c9d4SSatish Balay     goto done;
14319371c9d4SSatish Balay   }
14329566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX8", &flg));
14339371c9d4SSatish Balay   if (flg) {
14349371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
14359371c9d4SSatish Balay     goto done;
14369371c9d4SSatish Balay   }
14379566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg));
14389371c9d4SSatish Balay   if (flg) {
14399371c9d4SSatish Balay     *ct = DM_POLYTOPE_HEXAHEDRON;
14409371c9d4SSatish Balay     goto done;
14419371c9d4SSatish Balay   }
144298921bdaSJacob Faibussowitsch   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
14438f861fbcSMatthew G. Knepley done:
14448f861fbcSMatthew G. Knepley   PetscFunctionReturn(0);
14458f861fbcSMatthew G. Knepley }
14468f861fbcSMatthew G. Knepley #endif
14478f861fbcSMatthew G. Knepley 
1448552f7358SJed Brown /*@
1449a1cb98faSBarry Smith   DMPlexCreateExodus - Create a `DMPLEX` mesh from an ExodusII file ID.
1450552f7358SJed Brown 
1451d083f849SBarry Smith   Collective
1452552f7358SJed Brown 
1453552f7358SJed Brown   Input Parameters:
1454552f7358SJed Brown + comm  - The MPI communicator
1455552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1456552f7358SJed Brown - interpolate - Create faces and edges in the mesh
1457552f7358SJed Brown 
1458552f7358SJed Brown   Output Parameter:
1459a1cb98faSBarry Smith . dm  - The `DM` object representing the mesh
1460552f7358SJed Brown 
1461552f7358SJed Brown   Level: beginner
1462552f7358SJed Brown 
1463a1cb98faSBarry Smith .seealso: [](chapter_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `DMPLEX`, `DMCreate()`
1464552f7358SJed Brown @*/
1465d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
1466d71ae5a4SJacob Faibussowitsch {
1467552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1468552f7358SJed Brown   PetscMPIInt  num_proc, rank;
1469ae9eebabSLisandro Dalcin   DMLabel      cellSets = NULL, faceSets = NULL, vertSets = NULL;
1470552f7358SJed Brown   PetscSection coordSection;
1471552f7358SJed Brown   Vec          coordinates;
1472552f7358SJed Brown   PetscScalar *coords;
1473552f7358SJed Brown   PetscInt     coordSize, v;
1474552f7358SJed Brown   /* Read from ex_get_init() */
1475552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN + 1];
14765f80ce2aSJacob Faibussowitsch   int  dim = 0, dimEmbed = 0, numVertices = 0, numCells = 0;
1477552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
1478552f7358SJed Brown #endif
1479552f7358SJed Brown 
1480552f7358SJed Brown   PetscFunctionBegin;
1481552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
14829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
14839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
14849566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
14859566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1486a5b23f4aSJose E. Roman   /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */
1487dd400576SPatrick Sanan   if (rank == 0) {
14889566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(title, PETSC_MAX_PATH_LEN + 1));
1489792fecdfSBarry Smith     PetscCallExternal(ex_get_init, exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
149028b400f6SJacob Faibussowitsch     PetscCheck(num_cs, PETSC_COMM_SELF, PETSC_ERR_SUP, "Exodus file does not contain any cell set");
1491552f7358SJed Brown   }
14929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN + 1, MPI_CHAR, 0, comm));
14939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
14949566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)*dm, title));
14959566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, numCells + numVertices));
1496412e9a14SMatthew G. Knepley   /*   We do not want this label automatically computed, instead we compute it here */
14979566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "celltype"));
1498552f7358SJed Brown 
1499552f7358SJed Brown   /* Read cell sets information */
1500dd400576SPatrick Sanan   if (rank == 0) {
1501552f7358SJed Brown     PetscInt *cone;
1502e8f6893fSMatthew G. Knepley     int       c, cs, ncs, c_loc, v, v_loc;
1503552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
1504e8f6893fSMatthew G. Knepley     int *cs_id, *cs_order;
1505552f7358SJed Brown     /* Read from ex_get_elem_block() */
1506552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN + 1];
1507e8f6893fSMatthew G. Knepley     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1508552f7358SJed Brown     /* Read from ex_get_elem_conn() */
1509552f7358SJed Brown     int *cs_connect;
1510552f7358SJed Brown 
1511552f7358SJed Brown     /* Get cell sets IDs */
15129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order));
1513792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_ELEM_BLOCK, cs_id);
1514552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
1515552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1516e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
1517e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
15188f861fbcSMatthew G. Knepley       DMPolytopeType ct;
15198f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
15208f861fbcSMatthew G. Knepley 
15219566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1522792fecdfSBarry Smith       PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
15239566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
15248f861fbcSMatthew G. Knepley       dim = PetscMax(dim, DMPolytopeTypeGetDim(ct));
1525792fecdfSBarry 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);
15268f861fbcSMatthew G. Knepley       switch (ct) {
15278f861fbcSMatthew G. Knepley       case DM_POLYTOPE_TRI_PRISM:
1528e8f6893fSMatthew G. Knepley         cs_order[cs] = cs;
1529e8f6893fSMatthew G. Knepley         ++num_hybrid;
1530e8f6893fSMatthew G. Knepley         break;
1531e8f6893fSMatthew G. Knepley       default:
1532e8f6893fSMatthew G. Knepley         for (c = cs; c > cs - num_hybrid; --c) cs_order[c] = cs_order[c - 1];
1533e8f6893fSMatthew G. Knepley         cs_order[cs - num_hybrid] = cs;
1534e8f6893fSMatthew G. Knepley       }
1535e8f6893fSMatthew G. Knepley     }
1536552f7358SJed Brown     /* First set sizes */
1537e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
15388f861fbcSMatthew G. Knepley       DMPolytopeType ct;
15398f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
1540e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
15418f861fbcSMatthew G. Knepley 
15429566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1543792fecdfSBarry Smith       PetscCallExternal(ex_get_elem_type, exoid, cs_id[cs], elem_type);
15449566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
1545792fecdfSBarry 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);
1546552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
15479566063dSJacob Faibussowitsch         PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell));
15489566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCellType(*dm, c, ct));
1549552f7358SJed Brown       }
1550552f7358SJed Brown     }
15519566063dSJacob Faibussowitsch     for (v = numCells; v < numCells + numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
15529566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dm));
1553e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1554e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
1555792fecdfSBarry 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);
15569566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_vertex_per_cell * num_cell_in_set, &cs_connect, num_vertex_per_cell, &cone));
1557792fecdfSBarry Smith       PetscCallExternal(ex_get_conn, exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect, NULL, NULL);
1558eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1559552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1560636e64ffSStefano Zampini         DMPolytopeType ct;
1561636e64ffSStefano Zampini 
1562ad540459SPierre Jolivet         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) cone[v_loc] = cs_connect[v] + numCells - 1;
15639566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(*dm, c, &ct));
15649566063dSJacob Faibussowitsch         PetscCall(DMPlexInvertCell(ct, cone));
15659566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCone(*dm, c, cone));
15669566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]));
1567552f7358SJed Brown       }
15689566063dSJacob Faibussowitsch       PetscCall(PetscFree2(cs_connect, cone));
1569552f7358SJed Brown     }
15709566063dSJacob Faibussowitsch     PetscCall(PetscFree2(cs_id, cs_order));
1571552f7358SJed Brown   }
15728f861fbcSMatthew G. Knepley   {
15738f861fbcSMatthew G. Knepley     PetscInt ints[] = {dim, dimEmbed};
15748f861fbcSMatthew G. Knepley 
15759566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm));
15769566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dm, ints[0]));
15779566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(*dm, ints[1]));
15788f861fbcSMatthew G. Knepley     dim      = ints[0];
15798f861fbcSMatthew G. Knepley     dimEmbed = ints[1];
15808f861fbcSMatthew G. Knepley   }
15819566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
15829566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
1583552f7358SJed Brown   if (interpolate) {
15845fd9971aSMatthew G. Knepley     DM idm;
1585552f7358SJed Brown 
15869566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
15879566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1588552f7358SJed Brown     *dm = idm;
1589552f7358SJed Brown   }
1590552f7358SJed Brown 
1591552f7358SJed Brown   /* Create vertex set label */
1592dd400576SPatrick Sanan   if (rank == 0 && (num_vs > 0)) {
1593552f7358SJed Brown     int vs, v;
1594552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1595552f7358SJed Brown     int *vs_id;
1596552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1597f958083aSBlaise Bourdin     int num_vertex_in_set;
1598552f7358SJed Brown     /* Read from ex_get_node_set() */
1599552f7358SJed Brown     int *vs_vertex_list;
1600552f7358SJed Brown 
1601552f7358SJed Brown     /* Get vertex set ids */
16029566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_vs, &vs_id));
1603792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_NODE_SET, vs_id);
1604552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
1605792fecdfSBarry Smith       PetscCallExternal(ex_get_set_param, exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
16069566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list));
1607792fecdfSBarry Smith       PetscCallExternal(ex_get_set, exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);
160848a46eb9SPierre 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]));
16099566063dSJacob Faibussowitsch       PetscCall(PetscFree(vs_vertex_list));
1610552f7358SJed Brown     }
16119566063dSJacob Faibussowitsch     PetscCall(PetscFree(vs_id));
1612552f7358SJed Brown   }
1613552f7358SJed Brown   /* Read coordinates */
16149566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
16159566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
16169566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
16179566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
1618552f7358SJed Brown   for (v = numCells; v < numCells + numVertices; ++v) {
16199566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
16209566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
1621552f7358SJed Brown   }
16229566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
16239566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
16249566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
16259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject)coordinates, "coordinates"));
16269566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
16279566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, dimEmbed));
16289566063dSJacob Faibussowitsch   PetscCall(VecSetType(coordinates, VECSTANDARD));
16299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
1630dd400576SPatrick Sanan   if (rank == 0) {
1631e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1632552f7358SJed Brown 
16339566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVertices, &x, numVertices, &y, numVertices, &z));
1634792fecdfSBarry Smith     PetscCallExternal(ex_get_coord, exoid, x, y, z);
16358f861fbcSMatthew G. Knepley     if (dimEmbed > 0) {
16368f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 0] = x[v];
16370d644c17SKarl Rupp     }
16388f861fbcSMatthew G. Knepley     if (dimEmbed > 1) {
16398f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 1] = y[v];
16400d644c17SKarl Rupp     }
16418f861fbcSMatthew G. Knepley     if (dimEmbed > 2) {
16428f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v * dimEmbed + 2] = z[v];
16430d644c17SKarl Rupp     }
16449566063dSJacob Faibussowitsch     PetscCall(PetscFree3(x, y, z));
1645552f7358SJed Brown   }
16469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
16479566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
16489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
1649552f7358SJed Brown 
1650552f7358SJed Brown   /* Create side set label */
1651dd400576SPatrick Sanan   if (rank == 0 && interpolate && (num_fs > 0)) {
1652552f7358SJed Brown     int fs, f, voff;
1653552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1654552f7358SJed Brown     int *fs_id;
1655552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1656f958083aSBlaise Bourdin     int num_side_in_set;
1657552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1658552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
1659ef073283Smichael_afanasiev     /* Read side set labels */
1660ef073283Smichael_afanasiev     char   fs_name[MAX_STR_LENGTH + 1];
1661034d25ccSMatthew G. Knepley     size_t fs_name_len;
1662552f7358SJed Brown 
1663552f7358SJed Brown     /* Get side set ids */
16649566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_fs, &fs_id));
1665792fecdfSBarry Smith     PetscCallExternal(ex_get_ids, exoid, EX_SIDE_SET, fs_id);
1666552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
1667792fecdfSBarry Smith       PetscCallExternal(ex_get_set_param, exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);
16689566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_side_in_set, &fs_vertex_count_list, num_side_in_set * 4, &fs_vertex_list));
1669792fecdfSBarry Smith       PetscCallExternal(ex_get_side_set_node_list, exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);
1670ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1671ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1672034d25ccSMatthew G. Knepley       if (!fs_name_err) {
1673034d25ccSMatthew G. Knepley         PetscCall(PetscStrlen(fs_name, &fs_name_len));
1674034d25ccSMatthew G. Knepley         if (fs_name_len == 0) PetscCall(PetscStrncpy(fs_name, "Face Sets", MAX_STR_LENGTH));
1675034d25ccSMatthew G. Knepley       }
1676552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
16770298fd71SBarry Smith         const PetscInt *faces    = NULL;
1678552f7358SJed Brown         PetscInt        faceSize = fs_vertex_count_list[f], numFaces;
1679552f7358SJed Brown         PetscInt        faceVertices[4], v;
1680552f7358SJed Brown 
1681197f6770SJed Brown         PetscCheck(faceSize <= 4, comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %" PetscInt_FMT " > 4 vertices", faceSize);
1682ad540459SPierre Jolivet         for (v = 0; v < faceSize; ++v, ++voff) faceVertices[v] = fs_vertex_list[voff] + numCells - 1;
16839566063dSJacob Faibussowitsch         PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1684197f6770SJed Brown         PetscCheck(numFaces == 1, comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %" PetscInt_FMT " faces", f, fs, numFaces);
16859566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]));
1686ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1687034d25ccSMatthew G. Knepley         if (!fs_name_err) PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]));
16889566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1689552f7358SJed Brown       }
16909566063dSJacob Faibussowitsch       PetscCall(PetscFree2(fs_vertex_count_list, fs_vertex_list));
1691552f7358SJed Brown     }
16929566063dSJacob Faibussowitsch     PetscCall(PetscFree(fs_id));
1693552f7358SJed Brown   }
1694ae9eebabSLisandro Dalcin 
1695ae9eebabSLisandro Dalcin   { /* Create Cell/Face/Vertex Sets labels at all processes */
16969371c9d4SSatish Balay     enum {
16979371c9d4SSatish Balay       n = 3
16989371c9d4SSatish Balay     };
1699ae9eebabSLisandro Dalcin     PetscBool flag[n];
1700ae9eebabSLisandro Dalcin 
1701ae9eebabSLisandro Dalcin     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1702ae9eebabSLisandro Dalcin     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1703ae9eebabSLisandro Dalcin     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
17049566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
17059566063dSJacob Faibussowitsch     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
17069566063dSJacob Faibussowitsch     if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
17079566063dSJacob Faibussowitsch     if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
1708ae9eebabSLisandro Dalcin   }
1709b458e8f1SJose E. Roman   PetscFunctionReturn(0);
1710552f7358SJed Brown #else
1711552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1712552f7358SJed Brown #endif
1713552f7358SJed Brown }
1714