xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown 
4519f805aSKarl Rupp #if defined(PETSC_HAVE_EXODUSII)
5552f7358SJed Brown #include <netcdf.h>
6552f7358SJed Brown #include <exodusII.h>
7552f7358SJed Brown #endif
8552f7358SJed Brown 
91e50132fSMatthew G. Knepley #include <petsc/private/viewerimpl.h>
101e50132fSMatthew G. Knepley #include <petsc/private/viewerexodusiiimpl.h>
11e0266925SMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
1278569bb4SMatthew G. Knepley /*
131e50132fSMatthew G. Knepley   PETSC_VIEWER_EXODUSII_ - Creates an ExodusII PetscViewer shared by all processors in a communicator.
141e50132fSMatthew G. Knepley 
151e50132fSMatthew G. Knepley   Collective
161e50132fSMatthew G. Knepley 
171e50132fSMatthew G. Knepley   Input Parameter:
181e50132fSMatthew G. Knepley . comm - the MPI communicator to share the ExodusII PetscViewer
191e50132fSMatthew G. Knepley 
201e50132fSMatthew G. Knepley   Level: intermediate
211e50132fSMatthew G. Knepley 
221e50132fSMatthew G. Knepley   Notes:
231e50132fSMatthew G. Knepley     misses Fortran bindings
241e50132fSMatthew G. Knepley 
251e50132fSMatthew G. Knepley   Notes:
261e50132fSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_EXODUSII_ does not return
271e50132fSMatthew G. Knepley   an error code.  The GLVIS PetscViewer is usually used in the form
281e50132fSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_EXODUSII_(comm));
291e50132fSMatthew G. Knepley 
3000373969SVaclav Hapla .seealso: PetscViewerExodusIIOpen(), PetscViewerType, PetscViewerCreate(), PetscViewerDestroy()
311e50132fSMatthew G. Knepley */
321e50132fSMatthew G. Knepley PetscViewer PETSC_VIEWER_EXODUSII_(MPI_Comm comm)
331e50132fSMatthew G. Knepley {
341e50132fSMatthew G. Knepley   PetscViewer    viewer;
351e50132fSMatthew G. Knepley   PetscErrorCode ierr;
361e50132fSMatthew G. Knepley 
371e50132fSMatthew G. Knepley   PetscFunctionBegin;
381e50132fSMatthew G. Knepley   ierr = PetscViewerExodusIIOpen(comm, "mesh.exo", FILE_MODE_WRITE, &viewer);
3900373969SVaclav Hapla   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
401e50132fSMatthew G. Knepley   ierr = PetscObjectRegisterDestroy((PetscObject) viewer);
4100373969SVaclav Hapla   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_EXODUSII_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
421e50132fSMatthew G. Knepley   PetscFunctionReturn(viewer);
431e50132fSMatthew G. Knepley }
441e50132fSMatthew G. Knepley 
451e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerView_ExodusII(PetscViewer v, PetscViewer viewer)
461e50132fSMatthew G. Knepley {
476823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) v->data;
481e50132fSMatthew G. Knepley 
491e50132fSMatthew G. Knepley   PetscFunctionBegin;
509566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename:    %s\n", exo->filename));
519566063dSJacob Faibussowitsch   if (exo->exoid)    PetscCall(PetscViewerASCIIPrintf(viewer, "exoid:       %d\n", exo->exoid));
529566063dSJacob Faibussowitsch   if (exo->btype)    PetscCall(PetscViewerASCIIPrintf(viewer, "IO Mode:     %d\n", exo->btype));
539566063dSJacob Faibussowitsch   if (exo->order)    PetscCall(PetscViewerASCIIPrintf(viewer, "Mesh order:  %d\n", exo->order));
541e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
551e50132fSMatthew G. Knepley }
561e50132fSMatthew G. Knepley 
571e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v)
581e50132fSMatthew G. Knepley {
591e50132fSMatthew G. Knepley   PetscFunctionBegin;
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHead(PetscOptionsObject, "ExodusII PetscViewer Options"));
619566063dSJacob Faibussowitsch   PetscCall(PetscOptionsTail());
621e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
631e50132fSMatthew G. Knepley }
641e50132fSMatthew G. Knepley 
651e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer)
661e50132fSMatthew G. Knepley {
671e50132fSMatthew G. Knepley   PetscFunctionBegin;
681e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
691e50132fSMatthew G. Knepley }
701e50132fSMatthew G. Knepley 
711e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
721e50132fSMatthew G. Knepley {
731e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
741e50132fSMatthew G. Knepley 
751e50132fSMatthew G. Knepley   PetscFunctionBegin;
76a74df02fSJacob Faibussowitsch   if (exo->exoid >= 0) {PetscStackCallStandard(ex_close,exo->exoid);}
779566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo->filename));
789566063dSJacob Faibussowitsch   PetscCall(PetscFree(exo));
799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL));
809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL));
819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL));
829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetId_C",NULL));
839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetOrder_C",NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerSetOrder_C",NULL));
851e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
861e50132fSMatthew G. Knepley }
871e50132fSMatthew G. Knepley 
881e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
891e50132fSMatthew G. Knepley {
901e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
911e50132fSMatthew G. Knepley   PetscMPIInt           rank;
921e50132fSMatthew G. Knepley   int                   CPU_word_size, IO_word_size, EXO_mode;
936823f3c5SBlaise Bourdin   MPI_Info              mpi_info = MPI_INFO_NULL;
946823f3c5SBlaise Bourdin   float                 EXO_version;
951e50132fSMatthew G. Knepley 
969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank));
971e50132fSMatthew G. Knepley   CPU_word_size = sizeof(PetscReal);
981e50132fSMatthew G. Knepley   IO_word_size  = sizeof(PetscReal);
991e50132fSMatthew G. Knepley 
1001e50132fSMatthew G. Knepley   PetscFunctionBegin;
1016823f3c5SBlaise Bourdin   if (exo->exoid >= 0) {
102a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_close,exo->exoid);
1036823f3c5SBlaise Bourdin     exo->exoid = -1;
1046823f3c5SBlaise Bourdin   }
1059566063dSJacob Faibussowitsch   if (exo->filename) PetscCall(PetscFree(exo->filename));
1069566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &exo->filename));
1071e50132fSMatthew G. Knepley   switch (exo->btype) {
1081e50132fSMatthew G. Knepley   case FILE_MODE_READ:
1096823f3c5SBlaise Bourdin     EXO_mode = EX_READ;
1101e50132fSMatthew G. Knepley     break;
1111e50132fSMatthew G. Knepley   case FILE_MODE_APPEND:
1126823f3c5SBlaise Bourdin   case FILE_MODE_UPDATE:
1136823f3c5SBlaise Bourdin   case FILE_MODE_APPEND_UPDATE:
1146823f3c5SBlaise Bourdin     /* Will fail if the file does not already exist */
1156823f3c5SBlaise Bourdin     EXO_mode = EX_WRITE;
1161e50132fSMatthew G. Knepley     break;
1171e50132fSMatthew G. Knepley   case FILE_MODE_WRITE:
1186823f3c5SBlaise Bourdin     /*
1196823f3c5SBlaise Bourdin       exodus only allows writing geometry upon file creation, so we will let DMView create the file.
1206823f3c5SBlaise Bourdin     */
1216823f3c5SBlaise Bourdin     PetscFunctionReturn(0);
1221e50132fSMatthew G. Knepley     break;
1231e50132fSMatthew G. Knepley   default:
1241e50132fSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
1251e50132fSMatthew G. Knepley   }
1261e50132fSMatthew G. Knepley   #if defined(PETSC_USE_64BIT_INDICES)
1271e50132fSMatthew G. Knepley   EXO_mode += EX_ALL_INT64_API;
1281e50132fSMatthew G. Knepley   #endif
1296823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(name,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,PETSC_COMM_WORLD,mpi_info);
1302c71b3e2SJacob Faibussowitsch   PetscCheckFalse(exo->exoid < 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
1311e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1321e50132fSMatthew G. Knepley }
1331e50132fSMatthew G. Knepley 
1341e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name)
1351e50132fSMatthew G. Knepley {
1361e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1371e50132fSMatthew G. Knepley 
1381e50132fSMatthew G. Knepley   PetscFunctionBegin;
1391e50132fSMatthew G. Knepley   *name = exo->filename;
1401e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1411e50132fSMatthew G. Knepley }
1421e50132fSMatthew G. Knepley 
1431e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
1441e50132fSMatthew G. Knepley {
1451e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1461e50132fSMatthew G. Knepley 
1471e50132fSMatthew G. Knepley   PetscFunctionBegin;
1481e50132fSMatthew G. Knepley   exo->btype = type;
1491e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1501e50132fSMatthew G. Knepley }
1511e50132fSMatthew G. Knepley 
1521e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
1531e50132fSMatthew G. Knepley {
1541e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1551e50132fSMatthew G. Knepley 
1561e50132fSMatthew G. Knepley   PetscFunctionBegin;
1571e50132fSMatthew G. Knepley   *type = exo->btype;
1581e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1591e50132fSMatthew G. Knepley }
1601e50132fSMatthew G. Knepley 
1611e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid)
1621e50132fSMatthew G. Knepley {
1631e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1641e50132fSMatthew G. Knepley 
1651e50132fSMatthew G. Knepley   PetscFunctionBegin;
1661e50132fSMatthew G. Knepley   *exoid = exo->exoid;
1671e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1681e50132fSMatthew G. Knepley }
1691e50132fSMatthew G. Knepley 
1706823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order)
1711e50132fSMatthew G. Knepley {
1726823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1731e50132fSMatthew G. Knepley 
1741e50132fSMatthew G. Knepley   PetscFunctionBegin;
1756823f3c5SBlaise Bourdin   *order = exo->order;
1766823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
1776823f3c5SBlaise Bourdin }
1786823f3c5SBlaise Bourdin 
1796823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order)
1806823f3c5SBlaise Bourdin {
1816823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1826823f3c5SBlaise Bourdin 
1836823f3c5SBlaise Bourdin   PetscFunctionBegin;
1846823f3c5SBlaise Bourdin   exo->order = order;
1851e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1861e50132fSMatthew G. Knepley }
1871e50132fSMatthew G. Knepley 
1881e50132fSMatthew G. Knepley /*MC
18900373969SVaclav Hapla    PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
1901e50132fSMatthew G. Knepley 
19100373969SVaclav Hapla .seealso:  PetscViewerExodusIIOpen(), PetscViewerCreate(), PETSCVIEWERBINARY, PETSCVIEWERHDF5, DMView(),
1921e50132fSMatthew G. Knepley            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
1931e50132fSMatthew G. Knepley 
1941e50132fSMatthew G. Knepley   Level: beginner
1951e50132fSMatthew G. Knepley M*/
1961e50132fSMatthew G. Knepley 
1971e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
1981e50132fSMatthew G. Knepley {
1991e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo;
2001e50132fSMatthew G. Knepley 
2011e50132fSMatthew G. Knepley   PetscFunctionBegin;
2029566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(v,&exo));
2031e50132fSMatthew G. Knepley 
2041e50132fSMatthew G. Knepley   v->data                = (void*) exo;
2051e50132fSMatthew G. Knepley   v->ops->destroy        = PetscViewerDestroy_ExodusII;
2061e50132fSMatthew G. Knepley   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
2071e50132fSMatthew G. Knepley   v->ops->setup          = PetscViewerSetUp_ExodusII;
2081e50132fSMatthew G. Knepley   v->ops->view           = PetscViewerView_ExodusII;
2091e50132fSMatthew G. Knepley   v->ops->flush          = 0;
2101e50132fSMatthew G. Knepley   exo->btype             = (PetscFileMode) -1;
2111e50132fSMatthew G. Knepley   exo->filename          = 0;
2121e50132fSMatthew G. Knepley   exo->exoid             = -1;
2131e50132fSMatthew G. Knepley 
2149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII));
2159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII));
2169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII));
2179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII));
2189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetId_C",PetscViewerExodusIIGetId_ExodusII));
2199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerSetOrder_C",PetscViewerExodusIISetOrder_ExodusII));
2209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetOrder_C",PetscViewerExodusIIGetOrder_ExodusII));
2211e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
2221e50132fSMatthew G. Knepley }
2231e50132fSMatthew G. Knepley 
2241e50132fSMatthew G. Knepley /*
22578569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
22678569bb4SMatthew G. Knepley 
2276823f3c5SBlaise Bourdin   Collective
22878569bb4SMatthew G. Knepley 
22978569bb4SMatthew G. Knepley   Input Parameters:
23078569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
23178569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
23278569bb4SMatthew G. Knepley - name     - the name of the result
23378569bb4SMatthew G. Knepley 
23478569bb4SMatthew G. Knepley   Output Parameters:
23578569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
23678569bb4SMatthew G. Knepley 
23778569bb4SMatthew G. Knepley   Notes:
23878569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
23978569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
24078569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
24178569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
24278569bb4SMatthew G. Knepley 
24378569bb4SMatthew G. Knepley   Level: beginner
24478569bb4SMatthew G. Knepley 
245cd921712SStefano Zampini .seealso: EXOGetVarIndex(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO()
24678569bb4SMatthew G. Knepley */
2476823f3c5SBlaise Bourdin PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
24878569bb4SMatthew G. Knepley {
2496de834b4SMatthew G. Knepley   int            num_vars, i, j;
25078569bb4SMatthew G. Knepley   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
25178569bb4SMatthew G. Knepley   const int      num_suffix = 5;
25278569bb4SMatthew G. Knepley   char          *suffix[5];
25378569bb4SMatthew G. Knepley   PetscBool      flg;
25478569bb4SMatthew G. Knepley 
25578569bb4SMatthew G. Knepley   PetscFunctionBegin;
25678569bb4SMatthew G. Knepley   suffix[0] = (char *) "";
25778569bb4SMatthew G. Knepley   suffix[1] = (char *) "_X";
25878569bb4SMatthew G. Knepley   suffix[2] = (char *) "_XX";
25978569bb4SMatthew G. Knepley   suffix[3] = (char *) "_1";
26078569bb4SMatthew G. Knepley   suffix[4] = (char *) "_11";
26178569bb4SMatthew G. Knepley 
2626823f3c5SBlaise Bourdin   *varIndex = -1;
263a74df02fSJacob Faibussowitsch   PetscStackCallStandard(ex_get_variable_param,exoid, obj_type, &num_vars);
26478569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
265a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_variable_name,exoid, obj_type, i+1, var_name);
26678569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j) {
2679566063dSJacob Faibussowitsch       PetscCall(PetscStrncpy(ext_name, name, MAX_STR_LENGTH));
2689566063dSJacob Faibussowitsch       PetscCall(PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH));
2699566063dSJacob Faibussowitsch       PetscCall(PetscStrcasecmp(ext_name, var_name, &flg));
27078569bb4SMatthew G. Knepley       if (flg) {
27178569bb4SMatthew G. Knepley         *varIndex = i+1;
2726823f3c5SBlaise Bourdin       }
2736823f3c5SBlaise Bourdin     }
2746823f3c5SBlaise Bourdin   }
27578569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
27678569bb4SMatthew G. Knepley }
27778569bb4SMatthew G. Knepley 
27878569bb4SMatthew G. Knepley /*
2796823f3c5SBlaise Bourdin   DMView_PlexExodusII - Write a DM to disk in exodus format
28078569bb4SMatthew G. Knepley 
28178569bb4SMatthew G. Knepley   Collective on dm
28278569bb4SMatthew G. Knepley 
28378569bb4SMatthew G. Knepley   Input Parameters:
28478569bb4SMatthew G. Knepley + dm  - The dm to be written
2856823f3c5SBlaise Bourdin . viewer - an exodusII viewer
28678569bb4SMatthew G. Knepley 
28778569bb4SMatthew G. Knepley   Notes:
28878569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
28978569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
29078569bb4SMatthew G. Knepley 
2916823f3c5SBlaise Bourdin   If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices)
2926823f3c5SBlaise Bourdin   will be written.
29378569bb4SMatthew G. Knepley 
29478569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
29578569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
2966823f3c5SBlaise Bourdin   The order of the mesh shall be set using PetscViewerExodusIISetOrder
29778569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
29878569bb4SMatthew G. Knepley 
2996823f3c5SBlaise Bourdin   This function will only handle TRI, TET, QUAD, and HEX cells.
30078569bb4SMatthew G. Knepley   Level: beginner
30178569bb4SMatthew G. Knepley 
3026823f3c5SBlaise Bourdin .seealso:
30378569bb4SMatthew G. Knepley */
3046823f3c5SBlaise Bourdin PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer)
30578569bb4SMatthew G. Knepley {
30678569bb4SMatthew G. Knepley   enum ElemType {TRI, QUAD, TET, HEX};
30778569bb4SMatthew G. Knepley   MPI_Comm        comm;
3086823f3c5SBlaise Bourdin   PetscInt        degree; /* the order of the mesh */
30978569bb4SMatthew G. Knepley   /* Connectivity Variables */
31078569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
31178569bb4SMatthew G. Knepley   /* Cell Sets */
31278569bb4SMatthew G. Knepley   DMLabel         csLabel;
31378569bb4SMatthew G. Knepley   IS              csIS;
31478569bb4SMatthew G. Knepley   const PetscInt *csIdx;
31578569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
31678569bb4SMatthew G. Knepley   enum ElemType  *type;
317fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
31878569bb4SMatthew G. Knepley   /* Coordinate Variables */
31978569bb4SMatthew G. Knepley   DM              cdm;
3206823f3c5SBlaise Bourdin   PetscSection    coordSection;
32178569bb4SMatthew G. Knepley   Vec             coord;
32278569bb4SMatthew G. Knepley   PetscInt      **nodes;
323eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
32478569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
32578569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
32678569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
32778569bb4SMatthew G. Knepley   const char     *dmName;
328fe209ef9SBlaise Bourdin   PetscInt        nodesTriP1[4]  = {3,0,0,0};
329fe209ef9SBlaise Bourdin   PetscInt        nodesTriP2[4]  = {3,3,0,0};
330fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP1[4] = {4,0,0,0};
331fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP2[4] = {4,4,0,1};
332fe209ef9SBlaise Bourdin   PetscInt        nodesTetP1[4]  = {4,0,0,0};
333fe209ef9SBlaise Bourdin   PetscInt        nodesTetP2[4]  = {4,6,0,0};
334fe209ef9SBlaise Bourdin   PetscInt        nodesHexP1[4]  = {8,0,0,0};
335fe209ef9SBlaise Bourdin   PetscInt        nodesHexP2[4]  = {8,12,6,1};
33678569bb4SMatthew G. Knepley   PetscErrorCode  ierr;
33778569bb4SMatthew G. Knepley 
3386823f3c5SBlaise Bourdin   int             CPU_word_size, IO_word_size, EXO_mode;
3396823f3c5SBlaise Bourdin   float           EXO_version;
3406823f3c5SBlaise Bourdin 
3416823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
3426823f3c5SBlaise Bourdin 
34378569bb4SMatthew G. Knepley   PetscFunctionBegin;
3449566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
3459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
3469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
3476823f3c5SBlaise Bourdin 
3486823f3c5SBlaise Bourdin   /*
3496823f3c5SBlaise Bourdin     Creating coordSection is a collective operation so we do it somewhat out of sequence
3506823f3c5SBlaise Bourdin   */
3519566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &coordSection));
3529566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocalSetUp(dm));
3536823f3c5SBlaise Bourdin   if (!rank) {
3546823f3c5SBlaise Bourdin     switch (exo->btype) {
3556823f3c5SBlaise Bourdin     case FILE_MODE_READ:
3566823f3c5SBlaise Bourdin     case FILE_MODE_APPEND:
3576823f3c5SBlaise Bourdin     case FILE_MODE_UPDATE:
3586823f3c5SBlaise Bourdin     case FILE_MODE_APPEND_UPDATE:
3596823f3c5SBlaise Bourdin       /* exodusII does not allow writing geometry to an existing file */
36098921bdaSJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
3616823f3c5SBlaise Bourdin     case FILE_MODE_WRITE:
3626823f3c5SBlaise Bourdin       /* Create an empty file if one already exists*/
3636823f3c5SBlaise Bourdin       EXO_mode = EX_CLOBBER;
3646823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
3656823f3c5SBlaise Bourdin       EXO_mode += EX_ALL_INT64_API;
3666823f3c5SBlaise Bourdin #endif
3676823f3c5SBlaise Bourdin         CPU_word_size = sizeof(PetscReal);
3686823f3c5SBlaise Bourdin         IO_word_size  = sizeof(PetscReal);
3696823f3c5SBlaise Bourdin         exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
3702c71b3e2SJacob Faibussowitsch         PetscCheckFalse(exo->exoid < 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
3716823f3c5SBlaise Bourdin 
3726823f3c5SBlaise Bourdin       break;
3736823f3c5SBlaise Bourdin     default:
3746823f3c5SBlaise Bourdin       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3756823f3c5SBlaise Bourdin     }
3766823f3c5SBlaise Bourdin 
37778569bb4SMatthew G. Knepley     /* --- Get DM info --- */
3789566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject) dm, &dmName));
3799566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepth(dm, &depth));
3809566063dSJacob Faibussowitsch     PetscCall(DMGetDimension(dm, &dim));
3819566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
3829566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
3839566063dSJacob Faibussowitsch     PetscCall(DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd));
3849566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd));
3859566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
38678569bb4SMatthew G. Knepley     numCells    = cEnd - cStart;
38778569bb4SMatthew G. Knepley     numEdges    = eEnd - eStart;
38878569bb4SMatthew G. Knepley     numVertices = vEnd - vStart;
38978569bb4SMatthew G. Knepley     if (depth == 3) {numFaces = fEnd - fStart;}
39078569bb4SMatthew G. Knepley     else            {numFaces = 0;}
3919566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Cell Sets", &num_cs));
3929566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Vertex Sets", &num_vs));
3939566063dSJacob Faibussowitsch     PetscCall(DMGetLabelSize(dm, "Face Sets", &num_fs));
3949566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinatesLocal(dm, &coord));
3959566063dSJacob Faibussowitsch     PetscCall(DMGetCoordinateDM(dm, &cdm));
39678569bb4SMatthew G. Knepley     if (num_cs > 0) {
3979566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Cell Sets", &csLabel));
3989566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(csLabel, &csIS));
3999566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(csIS, &csIdx));
40078569bb4SMatthew G. Knepley     }
4019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &nodes));
40278569bb4SMatthew G. Knepley     /* Set element type for each block and compute total number of nodes */
4039566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_cs, &type));
40478569bb4SMatthew G. Knepley     numNodes = numVertices;
4056823f3c5SBlaise Bourdin 
4069566063dSJacob Faibussowitsch     PetscCall(PetscViewerExodusIIGetOrder(viewer, &degree));
40778569bb4SMatthew G. Knepley     if (degree == 2) {numNodes += numEdges;}
40878569bb4SMatthew G. Knepley     cellsNotInConnectivity = numCells;
40978569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
41078569bb4SMatthew G. Knepley       IS              stratumIS;
41178569bb4SMatthew G. Knepley       const PetscInt *cells;
41278569bb4SMatthew G. Knepley       PetscScalar    *xyz = NULL;
41378569bb4SMatthew G. Knepley       PetscInt        csSize, closureSize;
41478569bb4SMatthew G. Knepley 
4159566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4169566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4179566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
4189566063dSJacob Faibussowitsch       PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
41978569bb4SMatthew G. Knepley       switch (dim) {
42078569bb4SMatthew G. Knepley         case 2:
42178569bb4SMatthew G. Knepley           if      (closureSize == 3*dim) {type[cs] = TRI;}
42278569bb4SMatthew G. Knepley           else if (closureSize == 4*dim) {type[cs] = QUAD;}
42398921bdaSJacob Faibussowitsch           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
42478569bb4SMatthew G. Knepley           break;
42578569bb4SMatthew G. Knepley         case 3:
42678569bb4SMatthew G. Knepley           if      (closureSize == 4*dim) {type[cs] = TET;}
42778569bb4SMatthew G. Knepley           else if (closureSize == 8*dim) {type[cs] = HEX;}
42898921bdaSJacob Faibussowitsch           else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
42978569bb4SMatthew G. Knepley           break;
43098921bdaSJacob Faibussowitsch         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
43178569bb4SMatthew G. Knepley       }
43278569bb4SMatthew G. Knepley       if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
43378569bb4SMatthew G. Knepley       if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
4349566063dSJacob Faibussowitsch       PetscCall(DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz));
43578569bb4SMatthew G. Knepley       /* Set nodes and Element type */
43678569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
43778569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesTriP1;
43878569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTriP2;
43978569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
44078569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesQuadP1;
44178569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesQuadP2;
44278569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
44378569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesTetP1;
44478569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTetP2;
44578569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
44678569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesHexP1;
44778569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesHexP2;
44878569bb4SMatthew G. Knepley       }
44978569bb4SMatthew G. Knepley       /* Compute the number of cells not in the connectivity table */
45078569bb4SMatthew G. Knepley       cellsNotInConnectivity -= nodes[cs][3]*csSize;
45178569bb4SMatthew G. Knepley 
4529566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
4539566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
45478569bb4SMatthew G. Knepley     }
455a74df02fSJacob Faibussowitsch     if (num_cs > 0) {PetscStackCallStandard(ex_put_init,exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs);}
45678569bb4SMatthew G. Knepley     /* --- Connectivity --- */
45778569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
45878569bb4SMatthew G. Knepley       IS              stratumIS;
45978569bb4SMatthew G. Knepley       const PetscInt *cells;
460fe209ef9SBlaise Bourdin       PetscInt       *connect, off = 0;
461eae8a387SMatthew G. Knepley       PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
46278569bb4SMatthew G. Knepley       PetscInt        csSize, c, connectSize, closureSize;
46378569bb4SMatthew G. Knepley       char           *elem_type = NULL;
46478569bb4SMatthew G. Knepley       char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
46578569bb4SMatthew G. Knepley       char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
46678569bb4SMatthew G. Knepley       char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
46778569bb4SMatthew G. Knepley       char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
46878569bb4SMatthew G. Knepley 
4699566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
4709566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
4719566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
47278569bb4SMatthew G. Knepley       /* Set Element type */
47378569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
47478569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_tri3;
47578569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tri6;
47678569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
47778569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_quad4;
47878569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_quad9;
47978569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
48078569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_tet4;
48178569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tet10;
48278569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
48378569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_hex8;
48478569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_hex27;
48578569bb4SMatthew G. Knepley       }
48678569bb4SMatthew G. Knepley       connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
4879566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect));
488a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_block,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1);
48978569bb4SMatthew G. Knepley       /* Find number of vertices, edges, and faces in the closure */
49078569bb4SMatthew G. Knepley       verticesInClosure = nodes[cs][0];
49178569bb4SMatthew G. Knepley       if (depth > 1) {
49278569bb4SMatthew G. Knepley         if (dim == 2) {
4939566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &edgesInClosure));
49478569bb4SMatthew G. Knepley         } else if (dim == 3) {
49578569bb4SMatthew G. Knepley           PetscInt *closure = NULL;
49678569bb4SMatthew G. Knepley 
4979566063dSJacob Faibussowitsch           PetscCall(DMPlexGetConeSize(dm, cells[0], &facesInClosure));
4989566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
49978569bb4SMatthew G. Knepley           edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
5009566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure));
50178569bb4SMatthew G. Knepley         }
50278569bb4SMatthew G. Knepley       }
50378569bb4SMatthew G. Knepley       /* Get connectivity for each cell */
50478569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
50578569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
50678569bb4SMatthew G. Knepley         PetscInt  temp, i;
50778569bb4SMatthew G. Knepley 
5089566063dSJacob Faibussowitsch         PetscCall(DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
50978569bb4SMatthew G. Knepley         for (i = 0; i < connectSize; ++i) {
51078569bb4SMatthew G. Knepley           if (i < nodes[cs][0]) {/* Vertices */
511fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
512fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
51378569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
514fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
515fe209ef9SBlaise Bourdin             if (nodes[cs][2] == 0) connect[i+off] -= numFaces;
516fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
51778569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]) { /* Cells */
518fe209ef9SBlaise Bourdin             connect[i+off] = closure[0] + 1;
519fe209ef9SBlaise Bourdin             connect[i+off] -= skipCells;
52078569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]) { /* Faces */
521fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
522fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
52378569bb4SMatthew G. Knepley           } else {
524fe209ef9SBlaise Bourdin             connect[i+off] = -1;
52578569bb4SMatthew G. Knepley           }
52678569bb4SMatthew G. Knepley         }
52778569bb4SMatthew G. Knepley         /* Tetrahedra are inverted */
52878569bb4SMatthew G. Knepley         if (type[cs] == TET) {
529fe209ef9SBlaise Bourdin           temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp;
53078569bb4SMatthew G. Knepley           if (degree == 2) {
531e2c9586dSBlaise Bourdin             temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp;
532e2c9586dSBlaise Bourdin             temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp;
53378569bb4SMatthew G. Knepley           }
53478569bb4SMatthew G. Knepley         }
53578569bb4SMatthew G. Knepley         /* Hexahedra are inverted */
53678569bb4SMatthew G. Knepley         if (type[cs] == HEX) {
537fe209ef9SBlaise Bourdin           temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp;
53878569bb4SMatthew G. Knepley           if (degree == 2) {
539fe209ef9SBlaise Bourdin             temp = connect[8+off];  connect[8+off]  = connect[11+off]; connect[11+off] = temp;
540fe209ef9SBlaise Bourdin             temp = connect[9+off];  connect[9+off]  = connect[10+off]; connect[10+off] = temp;
541fe209ef9SBlaise Bourdin             temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp;
542fe209ef9SBlaise Bourdin             temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp;
54378569bb4SMatthew G. Knepley 
544fe209ef9SBlaise Bourdin             temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp;
545fe209ef9SBlaise Bourdin             temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp;
546fe209ef9SBlaise Bourdin             temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp;
547fe209ef9SBlaise Bourdin             temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp;
54878569bb4SMatthew G. Knepley 
549fe209ef9SBlaise Bourdin             temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp;
550fe209ef9SBlaise Bourdin             temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp;
551fe209ef9SBlaise Bourdin             temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp;
55278569bb4SMatthew G. Knepley           }
55378569bb4SMatthew G. Knepley         }
554fe209ef9SBlaise Bourdin         off += connectSize;
5559566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure));
55678569bb4SMatthew G. Knepley       }
557a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_conn,exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0);
55878569bb4SMatthew G. Knepley       skipCells += (nodes[cs][3] == 0)*csSize;
5599566063dSJacob Faibussowitsch       PetscCall(PetscFree(connect));
5609566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
5619566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
56278569bb4SMatthew G. Knepley     }
5639566063dSJacob Faibussowitsch     PetscCall(PetscFree(type));
56478569bb4SMatthew G. Knepley     /* --- Coordinates --- */
5659566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(coordSection, pStart, pEnd));
5661e50132fSMatthew G. Knepley     if (num_cs) {
56778569bb4SMatthew G. Knepley       for (d = 0; d < depth; ++d) {
5689566063dSJacob Faibussowitsch         PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
56978569bb4SMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
5709566063dSJacob Faibussowitsch           PetscCall(PetscSectionSetDof(coordSection, p, nodes[0][d] > 0));
57178569bb4SMatthew G. Knepley         }
57278569bb4SMatthew G. Knepley       }
5731e50132fSMatthew G. Knepley     }
57478569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
57578569bb4SMatthew G. Knepley       IS              stratumIS;
57678569bb4SMatthew G. Knepley       const PetscInt *cells;
57778569bb4SMatthew G. Knepley       PetscInt        csSize, c;
57878569bb4SMatthew G. Knepley 
5799566063dSJacob Faibussowitsch       PetscCall(DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS));
5809566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(stratumIS, &cells));
5819566063dSJacob Faibussowitsch       PetscCall(ISGetSize(stratumIS, &csSize));
58278569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
5839566063dSJacob Faibussowitsch         PetscCall(PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0));
58478569bb4SMatthew G. Knepley       }
5859566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(stratumIS, &cells));
5869566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&stratumIS));
58778569bb4SMatthew G. Knepley     }
58878569bb4SMatthew G. Knepley     if (num_cs > 0) {
5899566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(csIS, &csIdx));
5909566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&csIS));
59178569bb4SMatthew G. Knepley     }
5929566063dSJacob Faibussowitsch     PetscCall(PetscFree(nodes));
5939566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(coordSection));
59478569bb4SMatthew G. Knepley     if (numNodes > 0) {
59578569bb4SMatthew G. Knepley       const char  *coordNames[3] = {"x", "y", "z"};
596233c95e0SFrancesco Ballarin       PetscScalar *closure, *cval;
597233c95e0SFrancesco Ballarin       PetscReal   *coords;
59878569bb4SMatthew G. Knepley       PetscInt     hasDof, n = 0;
59978569bb4SMatthew G. Knepley 
6006823f3c5SBlaise Bourdin       /* There can't be more than 24 values in the closure of a point for the coord coordSection */
6019566063dSJacob Faibussowitsch       PetscCall(PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure));
6029566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinatesLocalNoncollective(dm, &coord));
6039566063dSJacob Faibussowitsch       PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
60478569bb4SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
6059566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(coordSection, p, &hasDof));
60678569bb4SMatthew G. Knepley         if (hasDof) {
60778569bb4SMatthew G. Knepley           PetscInt closureSize = 24, j;
60878569bb4SMatthew G. Knepley 
6099566063dSJacob Faibussowitsch           PetscCall(DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure));
61078569bb4SMatthew G. Knepley           for (d = 0; d < dim; ++d) {
61178569bb4SMatthew G. Knepley             cval[d] = 0.0;
61278569bb4SMatthew G. Knepley             for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
613233c95e0SFrancesco Ballarin             coords[d*numNodes+n] = PetscRealPart(cval[d]) * dim / closureSize;
61478569bb4SMatthew G. Knepley           }
61578569bb4SMatthew G. Knepley           ++n;
61678569bb4SMatthew G. Knepley         }
61778569bb4SMatthew G. Knepley       }
618a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_coord,exo->exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]);
6199566063dSJacob Faibussowitsch       PetscCall(PetscFree3(coords, cval, closure));
620a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_coord_names,exo->exoid, (char **) coordNames);
62178569bb4SMatthew G. Knepley     }
6226823f3c5SBlaise Bourdin 
623fe209ef9SBlaise Bourdin     /* --- Node Sets/Vertex Sets --- */
624fe209ef9SBlaise Bourdin     DMHasLabel(dm, "Vertex Sets", &hasLabel);
625fe209ef9SBlaise Bourdin     if (hasLabel) {
626fe209ef9SBlaise Bourdin       PetscInt        i, vs, vsSize;
627fe209ef9SBlaise Bourdin       const PetscInt *vsIdx, *vertices;
628fe209ef9SBlaise Bourdin       PetscInt       *nodeList;
629fe209ef9SBlaise Bourdin       IS              vsIS, stratumIS;
630fe209ef9SBlaise Bourdin       DMLabel         vsLabel;
6319566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Vertex Sets", &vsLabel));
6329566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(vsLabel, &vsIS));
6339566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(vsIS, &vsIdx));
634fe209ef9SBlaise Bourdin       for (vs=0; vs<num_vs; ++vs) {
6359566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS));
6369566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(stratumIS, &vertices));
6379566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &vsSize));
6389566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(vsSize, &nodeList));
639fe209ef9SBlaise Bourdin         for (i=0; i<vsSize; ++i) {
640fe209ef9SBlaise Bourdin           nodeList[i] = vertices[i] - skipCells + 1;
641fe209ef9SBlaise Bourdin         }
642a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0);
643a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_put_set,exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL);
6449566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(stratumIS, &vertices));
6459566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
6469566063dSJacob Faibussowitsch         PetscCall(PetscFree(nodeList));
647fe209ef9SBlaise Bourdin       }
6489566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(vsIS, &vsIdx));
6499566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&vsIS));
650fe209ef9SBlaise Bourdin     }
651fe209ef9SBlaise Bourdin     /* --- Side Sets/Face Sets --- */
6529566063dSJacob Faibussowitsch     PetscCall(DMHasLabel(dm, "Face Sets", &hasLabel));
653fe209ef9SBlaise Bourdin     if (hasLabel) {
654fe209ef9SBlaise Bourdin       PetscInt        i, j, fs, fsSize;
655fe209ef9SBlaise Bourdin       const PetscInt *fsIdx, *faces;
656fe209ef9SBlaise Bourdin       IS              fsIS, stratumIS;
657fe209ef9SBlaise Bourdin       DMLabel         fsLabel;
658fe209ef9SBlaise Bourdin       PetscInt        numPoints, *points;
659fe209ef9SBlaise Bourdin       PetscInt        elem_list_size = 0;
660fe209ef9SBlaise Bourdin       PetscInt       *elem_list, *elem_ind, *side_list;
661fe209ef9SBlaise Bourdin 
6629566063dSJacob Faibussowitsch       PetscCall(DMGetLabel(dm, "Face Sets", &fsLabel));
663fe209ef9SBlaise Bourdin       /* Compute size of Node List and Element List */
6649566063dSJacob Faibussowitsch       PetscCall(DMLabelGetValueIS(fsLabel, &fsIS));
6659566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(fsIS, &fsIdx));
666fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
6679566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
6689566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &fsSize));
669fe209ef9SBlaise Bourdin         elem_list_size += fsSize;
6709566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
671fe209ef9SBlaise Bourdin       }
672fe209ef9SBlaise Bourdin       ierr = PetscMalloc3(num_fs, &elem_ind,
673fe209ef9SBlaise Bourdin                           elem_list_size, &elem_list,
6749566063dSJacob Faibussowitsch                           elem_list_size, &side_list);PetscCall(ierr);
675fe209ef9SBlaise Bourdin       elem_ind[0] = 0;
676fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
6779566063dSJacob Faibussowitsch         PetscCall(DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS));
6789566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(stratumIS, &faces));
6799566063dSJacob Faibussowitsch         PetscCall(ISGetSize(stratumIS, &fsSize));
680fe209ef9SBlaise Bourdin         /* Set Parameters */
681a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_put_set_param,exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0);
682fe209ef9SBlaise Bourdin         /* Indices */
683fe209ef9SBlaise Bourdin         if (fs<num_fs-1) {
684fe209ef9SBlaise Bourdin           elem_ind[fs+1] = elem_ind[fs] + fsSize;
685fe209ef9SBlaise Bourdin         }
686fe209ef9SBlaise Bourdin 
687fe209ef9SBlaise Bourdin         for (i=0; i<fsSize; ++i) {
688fe209ef9SBlaise Bourdin           /* Element List */
689fe209ef9SBlaise Bourdin           points = NULL;
6909566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
691fe209ef9SBlaise Bourdin           elem_list[elem_ind[fs] + i] = points[2] +1;
6929566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points));
693fe209ef9SBlaise Bourdin 
694fe209ef9SBlaise Bourdin           /* Side List */
695fe209ef9SBlaise Bourdin           points = NULL;
6969566063dSJacob Faibussowitsch           PetscCall(DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points));
697fe209ef9SBlaise Bourdin           for (j=1; j<numPoints; ++j) {
698fe209ef9SBlaise Bourdin             if (points[j*2]==faces[i]) {break;}
699fe209ef9SBlaise Bourdin           }
700fe209ef9SBlaise Bourdin           /* Convert HEX sides */
701fe209ef9SBlaise Bourdin           if (numPoints == 27) {
702fe209ef9SBlaise Bourdin             if      (j == 1) {j = 5;}
703fe209ef9SBlaise Bourdin             else if (j == 2) {j = 6;}
704fe209ef9SBlaise Bourdin             else if (j == 3) {j = 1;}
705fe209ef9SBlaise Bourdin             else if (j == 4) {j = 3;}
706fe209ef9SBlaise Bourdin             else if (j == 5) {j = 2;}
707fe209ef9SBlaise Bourdin             else if (j == 6) {j = 4;}
708fe209ef9SBlaise Bourdin           }
709fe209ef9SBlaise Bourdin           /* Convert TET sides */
710fe209ef9SBlaise Bourdin           if (numPoints == 15) {
711fe209ef9SBlaise Bourdin             --j;
712fe209ef9SBlaise Bourdin             if (j == 0) {j = 4;}
713fe209ef9SBlaise Bourdin           }
714fe209ef9SBlaise Bourdin           side_list[elem_ind[fs] + i] = j;
7159566063dSJacob Faibussowitsch           PetscCall(DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points));
716fe209ef9SBlaise Bourdin 
717fe209ef9SBlaise Bourdin         }
7189566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(stratumIS, &faces));
7199566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&stratumIS));
720fe209ef9SBlaise Bourdin       }
7219566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(fsIS, &fsIdx));
7229566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&fsIS));
723fe209ef9SBlaise Bourdin 
724fe209ef9SBlaise Bourdin       /* Put side sets */
725fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
726a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_put_set,exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]);
727fe209ef9SBlaise Bourdin       }
7289566063dSJacob Faibussowitsch       PetscCall(PetscFree3(elem_ind,elem_list,side_list));
729fe209ef9SBlaise Bourdin     }
7306823f3c5SBlaise Bourdin     /*
7316823f3c5SBlaise Bourdin       close the exodus file
7326823f3c5SBlaise Bourdin     */
7336823f3c5SBlaise Bourdin     ex_close(exo->exoid);
7346823f3c5SBlaise Bourdin     exo->exoid = -1;
7356823f3c5SBlaise Bourdin   }
7369566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&coordSection));
7376823f3c5SBlaise Bourdin 
7386823f3c5SBlaise Bourdin   /*
7396823f3c5SBlaise Bourdin     reopen the file in parallel
7406823f3c5SBlaise Bourdin   */
7416823f3c5SBlaise Bourdin   EXO_mode = EX_WRITE;
7426823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
7436823f3c5SBlaise Bourdin   EXO_mode += EX_ALL_INT64_API;
7446823f3c5SBlaise Bourdin #endif
7456823f3c5SBlaise Bourdin   CPU_word_size = sizeof(PetscReal);
7466823f3c5SBlaise Bourdin   IO_word_size  = sizeof(PetscReal);
7476823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(exo->filename,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,comm,MPI_INFO_NULL);
7482c71b3e2SJacob Faibussowitsch   PetscCheckFalse(exo->exoid < 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
7496823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
7506823f3c5SBlaise Bourdin }
7516823f3c5SBlaise Bourdin 
7526823f3c5SBlaise Bourdin /*
7536823f3c5SBlaise Bourdin   VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
7546823f3c5SBlaise Bourdin 
7556823f3c5SBlaise Bourdin   Collective on v
7566823f3c5SBlaise Bourdin 
7576823f3c5SBlaise Bourdin   Input Parameters:
7586823f3c5SBlaise Bourdin + v  - The vector to be written
7596823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
7606823f3c5SBlaise Bourdin 
7616823f3c5SBlaise Bourdin   Notes:
7626823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
7636823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
7646823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
7656823f3c5SBlaise Bourdin   amongst all variables.
7666823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
7676823f3c5SBlaise Bourdin   possibly corrupting the file
7686823f3c5SBlaise Bourdin 
7696823f3c5SBlaise Bourdin   Level: beginner
7706823f3c5SBlaise Bourdin 
7716823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII(),VecView_PlexExodusII()
7726823f3c5SBlaise Bourdin @*/
7736823f3c5SBlaise Bourdin PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer)
7746823f3c5SBlaise Bourdin {
7756823f3c5SBlaise Bourdin   DM                 dm;
7766823f3c5SBlaise Bourdin   MPI_Comm           comm;
7776823f3c5SBlaise Bourdin   PetscMPIInt        rank;
7786823f3c5SBlaise Bourdin 
7796823f3c5SBlaise Bourdin   int                exoid,offsetN = 0, offsetZ = 0;
7806823f3c5SBlaise Bourdin   const char        *vecname;
7816823f3c5SBlaise Bourdin   PetscInt           step;
7826823f3c5SBlaise Bourdin 
7836823f3c5SBlaise Bourdin   PetscFunctionBegin;
7849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
7859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7869566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
7879566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
7889566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject) v, &vecname));
7896823f3c5SBlaise Bourdin 
7909566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL));
7919566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN));
7929566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ));
7932c71b3e2SJacob Faibussowitsch   PetscCheckFalse(offsetN <= 0 && offsetZ <= 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
7946823f3c5SBlaise Bourdin   if (offsetN > 0) {
7959566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN));
7966823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
7979566063dSJacob Faibussowitsch     PetscCall(VecViewPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ));
79898921bdaSJacob Faibussowitsch   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
7996823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
8006823f3c5SBlaise Bourdin }
8016823f3c5SBlaise Bourdin 
8026823f3c5SBlaise Bourdin /*
8036823f3c5SBlaise Bourdin   VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
8046823f3c5SBlaise Bourdin 
8056823f3c5SBlaise Bourdin   Collective on v
8066823f3c5SBlaise Bourdin 
8076823f3c5SBlaise Bourdin   Input Parameters:
8086823f3c5SBlaise Bourdin + v  - The vector to be written
8096823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
8106823f3c5SBlaise Bourdin 
8116823f3c5SBlaise Bourdin   Notes:
8126823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
8136823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
8146823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
8156823f3c5SBlaise Bourdin   amongst all variables.
8166823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
8176823f3c5SBlaise Bourdin   possibly corrupting the file
8186823f3c5SBlaise Bourdin 
8196823f3c5SBlaise Bourdin   Level: beginner
8206823f3c5SBlaise Bourdin 
8216823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII(),VecView_PlexExodusII()
8226823f3c5SBlaise Bourdin @*/
8236823f3c5SBlaise Bourdin PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer)
8246823f3c5SBlaise Bourdin {
8256823f3c5SBlaise Bourdin   DM                 dm;
8266823f3c5SBlaise Bourdin   MPI_Comm           comm;
8276823f3c5SBlaise Bourdin   PetscMPIInt        rank;
8286823f3c5SBlaise Bourdin 
8296823f3c5SBlaise Bourdin   int                exoid,offsetN = 0, offsetZ = 0;
8306823f3c5SBlaise Bourdin   const char        *vecname;
8316823f3c5SBlaise Bourdin   PetscInt           step;
8326823f3c5SBlaise Bourdin 
8336823f3c5SBlaise Bourdin   PetscFunctionBegin;
8349566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
8359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
8369566063dSJacob Faibussowitsch   PetscCall(PetscViewerExodusIIGetId(viewer,&exoid));
8379566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8389566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetName((PetscObject) v, &vecname));
8396823f3c5SBlaise Bourdin 
8409566063dSJacob Faibussowitsch   PetscCall(DMGetOutputSequenceNumber(dm,&step,NULL));
8419566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN));
8429566063dSJacob Faibussowitsch   PetscCall(EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ));
8432c71b3e2SJacob Faibussowitsch   PetscCheckFalse(offsetN <= 0 && offsetZ <= 0,comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. ", vecname);
8446823f3c5SBlaise Bourdin   if (offsetN > 0) {
8459566063dSJacob Faibussowitsch     PetscCall(VecLoadPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN));
8466823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
8479566063dSJacob Faibussowitsch     PetscCall(VecLoadPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ));
84898921bdaSJacob Faibussowitsch   } else SETERRQ(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. ", vecname);
84978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
85078569bb4SMatthew G. Knepley }
85178569bb4SMatthew G. Knepley 
85278569bb4SMatthew G. Knepley /*
85378569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
85478569bb4SMatthew G. Knepley 
85578569bb4SMatthew G. Knepley   Collective on v
85678569bb4SMatthew G. Knepley 
85778569bb4SMatthew G. Knepley   Input Parameters:
85878569bb4SMatthew G. Knepley + v  - The vector to be written
85978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
8606823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
8616823f3c5SBlaise Bourdin - offset - the location of the variable in the file
86278569bb4SMatthew G. Knepley 
86378569bb4SMatthew G. Knepley   Notes:
86478569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
86578569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
86678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
86778569bb4SMatthew G. Knepley   amongst all nodal variables.
86878569bb4SMatthew G. Knepley 
86978569bb4SMatthew G. Knepley   Level: beginner
87078569bb4SMatthew G. Knepley 
8716823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
87278569bb4SMatthew G. Knepley @*/
8736823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
87478569bb4SMatthew G. Knepley {
87578569bb4SMatthew G. Knepley   MPI_Comm           comm;
87678569bb4SMatthew G. Knepley   PetscMPIInt        size;
87778569bb4SMatthew G. Knepley   DM                 dm;
87878569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
87922a7544eSVaclav Hapla   const PetscScalar *varray;
88078569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
88178569bb4SMatthew G. Knepley   PetscBool          useNatural;
88278569bb4SMatthew G. Knepley 
88378569bb4SMatthew G. Knepley   PetscFunctionBegin;
8849566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
8859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
8869566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
8879566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
88878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
88978569bb4SMatthew G. Knepley   if (useNatural) {
8909566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &vNatural));
8919566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
8929566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
89378569bb4SMatthew G. Knepley   } else {
89478569bb4SMatthew G. Knepley     vNatural = v;
89578569bb4SMatthew G. Knepley   }
8966823f3c5SBlaise Bourdin 
89778569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
89878569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
89978569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
9009566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
9019566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
90278569bb4SMatthew G. Knepley   if (bs == 1) {
9039566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(vNatural, &varray));
904a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);
9059566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(vNatural, &varray));
90678569bb4SMatthew G. Knepley   } else {
90778569bb4SMatthew G. Knepley     IS       compIS;
90878569bb4SMatthew G. Knepley     PetscInt c;
90978569bb4SMatthew G. Knepley 
9109566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
91178569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
9129566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
9139566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
9149566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vComp, &varray));
915a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);
9169566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vComp, &varray));
9179566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
91878569bb4SMatthew G. Knepley     }
9199566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
92078569bb4SMatthew G. Knepley   }
9219566063dSJacob Faibussowitsch   if (useNatural) PetscCall(DMRestoreGlobalVector(dm, &vNatural));
92278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
92378569bb4SMatthew G. Knepley }
92478569bb4SMatthew G. Knepley 
92578569bb4SMatthew G. Knepley /*
92678569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
92778569bb4SMatthew G. Knepley 
92878569bb4SMatthew G. Knepley   Collective on v
92978569bb4SMatthew G. Knepley 
93078569bb4SMatthew G. Knepley   Input Parameters:
93178569bb4SMatthew G. Knepley + v  - The vector to be written
93278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
9336823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
9346823f3c5SBlaise Bourdin - offset - the location of the variable in the file
93578569bb4SMatthew G. Knepley 
93678569bb4SMatthew G. Knepley   Notes:
93778569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
93878569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
93978569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
94078569bb4SMatthew G. Knepley   amongst all nodal variables.
94178569bb4SMatthew G. Knepley 
94278569bb4SMatthew G. Knepley   Level: beginner
94378569bb4SMatthew G. Knepley 
9446823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
94578569bb4SMatthew G. Knepley */
9466823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
94778569bb4SMatthew G. Knepley {
94878569bb4SMatthew G. Knepley   MPI_Comm       comm;
94978569bb4SMatthew G. Knepley   PetscMPIInt    size;
95078569bb4SMatthew G. Knepley   DM             dm;
95178569bb4SMatthew G. Knepley   Vec            vNatural, vComp;
95222a7544eSVaclav Hapla   PetscScalar   *varray;
95378569bb4SMatthew G. Knepley   PetscInt       xs, xe, bs;
95478569bb4SMatthew G. Knepley   PetscBool      useNatural;
95578569bb4SMatthew G. Knepley 
95678569bb4SMatthew G. Knepley   PetscFunctionBegin;
9579566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
9589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
9599566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v,&dm));
9609566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
96178569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
9629566063dSJacob Faibussowitsch   if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural));
96378569bb4SMatthew G. Knepley   else            {vNatural = v;}
9646823f3c5SBlaise Bourdin 
96578569bb4SMatthew G. Knepley   /* Read local chunk from the file */
9669566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
9679566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
96878569bb4SMatthew G. Knepley   if (bs == 1) {
9699566063dSJacob Faibussowitsch     PetscCall(VecGetArray(vNatural, &varray));
970a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray);
9719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vNatural, &varray));
97278569bb4SMatthew G. Knepley   } else {
97378569bb4SMatthew G. Knepley     IS       compIS;
97478569bb4SMatthew G. Knepley     PetscInt c;
97578569bb4SMatthew G. Knepley 
9769566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
97778569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
9789566063dSJacob Faibussowitsch       PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
9799566063dSJacob Faibussowitsch       PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
9809566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vComp, &varray));
981a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray);
9829566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vComp, &varray));
9839566063dSJacob Faibussowitsch       PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
98478569bb4SMatthew G. Knepley     }
9859566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&compIS));
98678569bb4SMatthew G. Knepley   }
98778569bb4SMatthew G. Knepley   if (useNatural) {
9889566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
9899566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
9909566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &vNatural));
99178569bb4SMatthew G. Knepley   }
99278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
99378569bb4SMatthew G. Knepley }
99478569bb4SMatthew G. Knepley 
99578569bb4SMatthew G. Knepley /*
99678569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
99778569bb4SMatthew G. Knepley 
99878569bb4SMatthew G. Knepley   Collective on v
99978569bb4SMatthew G. Knepley 
100078569bb4SMatthew G. Knepley   Input Parameters:
100178569bb4SMatthew G. Knepley + v  - The vector to be written
100278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
10036823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
10046823f3c5SBlaise Bourdin - offset - the location of the variable in the file
100578569bb4SMatthew G. Knepley 
100678569bb4SMatthew G. Knepley   Notes:
100778569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
100878569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
100978569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
101078569bb4SMatthew G. Knepley   amongst all zonal variables.
101178569bb4SMatthew G. Knepley 
101278569bb4SMatthew G. Knepley   Level: beginner
101378569bb4SMatthew G. Knepley 
10146823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
101578569bb4SMatthew G. Knepley */
10166823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
101778569bb4SMatthew G. Knepley {
101878569bb4SMatthew G. Knepley   MPI_Comm          comm;
101978569bb4SMatthew G. Knepley   PetscMPIInt       size;
102078569bb4SMatthew G. Knepley   DM                dm;
102178569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
102222a7544eSVaclav Hapla   const PetscScalar *varray;
102378569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
102478569bb4SMatthew G. Knepley   PetscBool         useNatural;
102578569bb4SMatthew G. Knepley   IS                compIS;
102678569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
102778569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
102878569bb4SMatthew G. Knepley 
102978569bb4SMatthew G. Knepley   PetscFunctionBegin;
10309566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
10319566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
10329566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
10339566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
103478569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
103578569bb4SMatthew G. Knepley   if (useNatural) {
10369566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &vNatural));
10379566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalBegin(dm, v, vNatural));
10389566063dSJacob Faibussowitsch     PetscCall(DMPlexGlobalToNaturalEnd(dm, v, vNatural));
103978569bb4SMatthew G. Knepley   } else {
104078569bb4SMatthew G. Knepley     vNatural = v;
104178569bb4SMatthew G. Knepley   }
10426823f3c5SBlaise Bourdin 
104378569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
104478569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
104578569bb4SMatthew G. Knepley      We assume that they are stored sequentially
104678569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1047a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
104878569bb4SMatthew G. Knepley      to figure out what to save where. */
104978569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
10509566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1051a74df02fSJacob Faibussowitsch   PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID);
105278569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
105378569bb4SMatthew G. Knepley     ex_block block;
105478569bb4SMatthew G. Knepley 
105578569bb4SMatthew G. Knepley     block.id   = csID[set];
105678569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1057a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_block_param,exoid, &block);
105878569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
105978569bb4SMatthew G. Knepley   }
10609566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
10619566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
10629566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
106378569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
106478569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
106578569bb4SMatthew G. Knepley 
106678569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
106778569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
106878569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
106978569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
107078569bb4SMatthew G. Knepley     if (bs == 1) {
10719566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(vNatural, &varray));
1072a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);
10739566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(vNatural, &varray));
107478569bb4SMatthew G. Knepley     } else {
107578569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
10769566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
10779566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
10789566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(vComp, &varray));
1079a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_put_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]);
10809566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(vComp, &varray));
10819566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
108278569bb4SMatthew G. Knepley       }
108378569bb4SMatthew G. Knepley     }
108478569bb4SMatthew G. Knepley     csxs += csSize[set];
108578569bb4SMatthew G. Knepley   }
10869566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
10879566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
10889566063dSJacob Faibussowitsch   if (useNatural) PetscCall(DMRestoreGlobalVector(dm,&vNatural));
108978569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
109078569bb4SMatthew G. Knepley }
109178569bb4SMatthew G. Knepley 
109278569bb4SMatthew G. Knepley /*
109378569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
109478569bb4SMatthew G. Knepley 
109578569bb4SMatthew G. Knepley   Collective on v
109678569bb4SMatthew G. Knepley 
109778569bb4SMatthew G. Knepley   Input Parameters:
109878569bb4SMatthew G. Knepley + v  - The vector to be written
109978569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
11006823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
11016823f3c5SBlaise Bourdin - offset - the location of the variable in the file
110278569bb4SMatthew G. Knepley 
110378569bb4SMatthew G. Knepley   Notes:
110478569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
110578569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
110678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
110778569bb4SMatthew G. Knepley   amongst all zonal variables.
110878569bb4SMatthew G. Knepley 
110978569bb4SMatthew G. Knepley   Level: beginner
111078569bb4SMatthew G. Knepley 
11116823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
111278569bb4SMatthew G. Knepley */
11136823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
111478569bb4SMatthew G. Knepley {
111578569bb4SMatthew G. Knepley   MPI_Comm          comm;
111678569bb4SMatthew G. Knepley   PetscMPIInt       size;
111778569bb4SMatthew G. Knepley   DM                dm;
111878569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
111922a7544eSVaclav Hapla   PetscScalar      *varray;
112078569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
112178569bb4SMatthew G. Knepley   PetscBool         useNatural;
112278569bb4SMatthew G. Knepley   IS                compIS;
112378569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
112478569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
112578569bb4SMatthew G. Knepley 
112678569bb4SMatthew G. Knepley   PetscFunctionBegin;
11279566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v,&comm));
11289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
11299566063dSJacob Faibussowitsch   PetscCall(VecGetDM(v, &dm));
11309566063dSJacob Faibussowitsch   PetscCall(DMGetUseNatural(dm, &useNatural));
113178569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
11329566063dSJacob Faibussowitsch   if (useNatural) PetscCall(DMGetGlobalVector(dm,&vNatural));
113378569bb4SMatthew G. Knepley   else            {vNatural = v;}
11346823f3c5SBlaise Bourdin 
113578569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
113678569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
113778569bb4SMatthew G. Knepley      We assume that they are stored sequentially
113878569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
1139a5b23f4aSJose E. Roman      but once the vector has been reordered to natural size, we cannot use the label information
114078569bb4SMatthew G. Knepley      to figure out what to save where. */
114178569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
11429566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(numCS, &csID, numCS, &csSize));
1143a74df02fSJacob Faibussowitsch   PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, csID);
114478569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
114578569bb4SMatthew G. Knepley     ex_block block;
114678569bb4SMatthew G. Knepley 
114778569bb4SMatthew G. Knepley     block.id   = csID[set];
114878569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
1149a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_block_param,exoid, &block);
115078569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
115178569bb4SMatthew G. Knepley   }
11529566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(vNatural, &xs, &xe));
11539566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(vNatural, &bs));
11549566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS));
115578569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
115678569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
115778569bb4SMatthew G. Knepley 
115878569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
115978569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
116078569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
116178569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
116278569bb4SMatthew G. Knepley     if (bs == 1) {
11639566063dSJacob Faibussowitsch       PetscCall(VecGetArray(vNatural, &varray));
1164a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]);
11659566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(vNatural, &varray));
116678569bb4SMatthew G. Knepley     } else {
116778569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
11689566063dSJacob Faibussowitsch         PetscCall(ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs));
11699566063dSJacob Faibussowitsch         PetscCall(VecGetSubVector(vNatural, compIS, &vComp));
11709566063dSJacob Faibussowitsch         PetscCall(VecGetArray(vComp, &varray));
1171a74df02fSJacob Faibussowitsch         PetscStackCallStandard(ex_get_partial_var,exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]);
11729566063dSJacob Faibussowitsch         PetscCall(VecRestoreArray(vComp, &varray));
11739566063dSJacob Faibussowitsch         PetscCall(VecRestoreSubVector(vNatural, compIS, &vComp));
117478569bb4SMatthew G. Knepley       }
117578569bb4SMatthew G. Knepley     }
117678569bb4SMatthew G. Knepley     csxs += csSize[set];
117778569bb4SMatthew G. Knepley   }
11789566063dSJacob Faibussowitsch   PetscCall(PetscFree2(csID, csSize));
11799566063dSJacob Faibussowitsch   if (bs > 1) PetscCall(ISDestroy(&compIS));
118078569bb4SMatthew G. Knepley   if (useNatural) {
11819566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalBegin(dm, vNatural, v));
11829566063dSJacob Faibussowitsch     PetscCall(DMPlexNaturalToGlobalEnd(dm, vNatural, v));
11839566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &vNatural));
118478569bb4SMatthew G. Knepley   }
118578569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
118678569bb4SMatthew G. Knepley }
1187b53c8456SSatish Balay #endif
118878569bb4SMatthew G. Knepley 
11891e50132fSMatthew G. Knepley /*@
11901e50132fSMatthew G. Knepley   PetscViewerExodusIIGetId - Get the file id of the ExodusII file
11911e50132fSMatthew G. Knepley 
11921e50132fSMatthew G. Knepley   Logically Collective on PetscViewer
11931e50132fSMatthew G. Knepley 
11941e50132fSMatthew G. Knepley   Input Parameter:
11951e50132fSMatthew G. Knepley .  viewer - the PetscViewer
11961e50132fSMatthew G. Knepley 
11971e50132fSMatthew G. Knepley   Output Parameter:
1198d8d19677SJose E. Roman .  exoid - The ExodusII file id
11991e50132fSMatthew G. Knepley 
12001e50132fSMatthew G. Knepley   Level: intermediate
12011e50132fSMatthew G. Knepley 
12021e50132fSMatthew G. Knepley .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
12031e50132fSMatthew G. Knepley @*/
12041e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid)
12051e50132fSMatthew G. Knepley {
12061e50132fSMatthew G. Knepley   PetscFunctionBegin;
12071e50132fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1208*cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetId_C",(PetscViewer,int*),(viewer,exoid));
12096823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12106823f3c5SBlaise Bourdin }
12116823f3c5SBlaise Bourdin 
12126823f3c5SBlaise Bourdin /*@
12136823f3c5SBlaise Bourdin    PetscViewerExodusIISetOrder - Set the elements order in the exodusII file.
12146823f3c5SBlaise Bourdin 
12156823f3c5SBlaise Bourdin    Collective
12166823f3c5SBlaise Bourdin 
12176823f3c5SBlaise Bourdin    Input Parameters:
12186823f3c5SBlaise Bourdin +  viewer - the viewer
121998a6a78aSPierre Jolivet -  order - elements order
12206823f3c5SBlaise Bourdin 
12216823f3c5SBlaise Bourdin    Output Parameter:
12226823f3c5SBlaise Bourdin 
12236823f3c5SBlaise Bourdin    Level: beginner
12246823f3c5SBlaise Bourdin 
12256823f3c5SBlaise Bourdin    Note:
12266823f3c5SBlaise Bourdin 
12276823f3c5SBlaise Bourdin .seealso: PetscViewerExodusIIGetId(), PetscViewerExodusIIGetOrder(), PetscViewerExodusIISetOrder()
12286823f3c5SBlaise Bourdin @*/
12296823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order)
12306823f3c5SBlaise Bourdin {
12316823f3c5SBlaise Bourdin   PetscFunctionBegin;
12326823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1233*cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerSetOrder_C",(PetscViewer,PetscInt),(viewer,order));
12346823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12356823f3c5SBlaise Bourdin }
12366823f3c5SBlaise Bourdin 
12376823f3c5SBlaise Bourdin /*@
12386823f3c5SBlaise Bourdin    PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file.
12396823f3c5SBlaise Bourdin 
12406823f3c5SBlaise Bourdin    Collective
12416823f3c5SBlaise Bourdin 
12426823f3c5SBlaise Bourdin    Input Parameters:
12436823f3c5SBlaise Bourdin +  viewer - the viewer
124498a6a78aSPierre Jolivet -  order - elements order
12456823f3c5SBlaise Bourdin 
12466823f3c5SBlaise Bourdin    Output Parameter:
12476823f3c5SBlaise Bourdin 
12486823f3c5SBlaise Bourdin    Level: beginner
12496823f3c5SBlaise Bourdin 
12506823f3c5SBlaise Bourdin    Note:
12516823f3c5SBlaise Bourdin 
12526823f3c5SBlaise Bourdin .seealso: PetscViewerExodusIIGetId(), PetscViewerExodusIIGetOrder(), PetscViewerExodusIISetOrder()
12536823f3c5SBlaise Bourdin @*/
12546823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order)
12556823f3c5SBlaise Bourdin {
12566823f3c5SBlaise Bourdin   PetscFunctionBegin;
12576823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1258*cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerGetOrder_C",(PetscViewer,PetscInt*),(viewer,order));
12596823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12606823f3c5SBlaise Bourdin }
12616823f3c5SBlaise Bourdin 
12626823f3c5SBlaise Bourdin /*@C
12636823f3c5SBlaise Bourdin    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
12646823f3c5SBlaise Bourdin 
12656823f3c5SBlaise Bourdin    Collective
12666823f3c5SBlaise Bourdin 
12676823f3c5SBlaise Bourdin    Input Parameters:
12686823f3c5SBlaise Bourdin +  comm - MPI communicator
12696823f3c5SBlaise Bourdin .  name - name of file
12706823f3c5SBlaise Bourdin -  type - type of file
12716823f3c5SBlaise Bourdin $    FILE_MODE_WRITE - create new file for binary output
12726823f3c5SBlaise Bourdin $    FILE_MODE_READ - open existing file for binary input
12736823f3c5SBlaise Bourdin $    FILE_MODE_APPEND - open existing file for binary output
12746823f3c5SBlaise Bourdin 
12756823f3c5SBlaise Bourdin    Output Parameter:
12766823f3c5SBlaise Bourdin .  exo - PetscViewer for Exodus II input/output to use with the specified file
12776823f3c5SBlaise Bourdin 
12786823f3c5SBlaise Bourdin    Level: beginner
12796823f3c5SBlaise Bourdin 
12806823f3c5SBlaise Bourdin    Note:
12816823f3c5SBlaise Bourdin    This PetscViewer should be destroyed with PetscViewerDestroy().
12826823f3c5SBlaise Bourdin 
12836823f3c5SBlaise Bourdin .seealso: PetscViewerPushFormat(), PetscViewerDestroy(),
12846823f3c5SBlaise Bourdin           DMLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
12856823f3c5SBlaise Bourdin @*/
12866823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo)
12876823f3c5SBlaise Bourdin {
12886823f3c5SBlaise Bourdin   PetscFunctionBegin;
12899566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, exo));
12909566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*exo, PETSCVIEWEREXODUSII));
12919566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*exo, type));
12929566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*exo, name));
12939566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*exo));
12941e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
12951e50132fSMatthew G. Knepley }
12961e50132fSMatthew G. Knepley 
129733751fbdSMichael Lange /*@C
1298eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
129933751fbdSMichael Lange 
1300d083f849SBarry Smith   Collective
130133751fbdSMichael Lange 
130233751fbdSMichael Lange   Input Parameters:
130333751fbdSMichael Lange + comm  - The MPI communicator
130433751fbdSMichael Lange . filename - The name of the ExodusII file
130533751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
130633751fbdSMichael Lange 
130733751fbdSMichael Lange   Output Parameter:
130833751fbdSMichael Lange . dm  - The DM object representing the mesh
130933751fbdSMichael Lange 
131033751fbdSMichael Lange   Level: beginner
131133751fbdSMichael Lange 
131233751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
131333751fbdSMichael Lange @*/
131433751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
131533751fbdSMichael Lange {
131633751fbdSMichael Lange   PetscMPIInt    rank;
131733751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1318e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
131933751fbdSMichael Lange   float version;
132033751fbdSMichael Lange #endif
132133751fbdSMichael Lange 
132233751fbdSMichael Lange   PetscFunctionBegin;
132333751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
13249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
132533751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1326dd400576SPatrick Sanan   if (rank == 0) {
132733751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
13282c71b3e2SJacob Faibussowitsch     PetscCheckFalse(exoid <= 0,PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
132933751fbdSMichael Lange   }
13309566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateExodus(comm, exoid, interpolate, dm));
1331a74df02fSJacob Faibussowitsch   if (rank == 0) {PetscStackCallStandard(ex_close,exoid);}
1332b458e8f1SJose E. Roman   PetscFunctionReturn(0);
133333751fbdSMichael Lange #else
133433751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
133533751fbdSMichael Lange #endif
133633751fbdSMichael Lange }
133733751fbdSMichael Lange 
13388f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
13396823f3c5SBlaise Bourdin static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct)
13408f861fbcSMatthew G. Knepley {
13418f861fbcSMatthew G. Knepley   PetscBool      flg;
13428f861fbcSMatthew G. Knepley 
13438f861fbcSMatthew G. Knepley   PetscFunctionBegin;
13448f861fbcSMatthew G. Knepley   *ct = DM_POLYTOPE_UNKNOWN;
13459566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI", &flg));
13468f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
13479566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TRI3", &flg));
13488f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
13499566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD", &flg));
13508f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13519566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "QUAD4", &flg));
13528f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13539566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "SHELL4", &flg));
13548f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13559566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TETRA", &flg));
13568f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
13579566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "TET4", &flg));
13588f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
13599566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "WEDGE", &flg));
13608f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRI_PRISM; goto done;}
13619566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX", &flg));
13628f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
13639566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEX8", &flg));
13648f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
13659566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(elem_type, "HEXAHEDRON", &flg));
13668f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
136798921bdaSJacob Faibussowitsch   SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
13688f861fbcSMatthew G. Knepley   done:
13698f861fbcSMatthew G. Knepley   PetscFunctionReturn(0);
13708f861fbcSMatthew G. Knepley }
13718f861fbcSMatthew G. Knepley #endif
13728f861fbcSMatthew G. Knepley 
1373552f7358SJed Brown /*@
137433751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
1375552f7358SJed Brown 
1376d083f849SBarry Smith   Collective
1377552f7358SJed Brown 
1378552f7358SJed Brown   Input Parameters:
1379552f7358SJed Brown + comm  - The MPI communicator
1380552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1381552f7358SJed Brown - interpolate - Create faces and edges in the mesh
1382552f7358SJed Brown 
1383552f7358SJed Brown   Output Parameter:
1384552f7358SJed Brown . dm  - The DM object representing the mesh
1385552f7358SJed Brown 
1386552f7358SJed Brown   Level: beginner
1387552f7358SJed Brown 
1388e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
1389552f7358SJed Brown @*/
1390552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
1391552f7358SJed Brown {
1392552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1393552f7358SJed Brown   PetscMPIInt    num_proc, rank;
1394ae9eebabSLisandro Dalcin   DMLabel        cellSets = NULL, faceSets = NULL, vertSets = NULL;
1395552f7358SJed Brown   PetscSection   coordSection;
1396552f7358SJed Brown   Vec            coordinates;
1397552f7358SJed Brown   PetscScalar    *coords;
1398552f7358SJed Brown   PetscInt       coordSize, v;
1399552f7358SJed Brown   /* Read from ex_get_init() */
1400552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
14015f80ce2aSJacob Faibussowitsch   int  dim    = 0, dimEmbed = 0, numVertices = 0, numCells = 0;
1402552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
1403552f7358SJed Brown #endif
1404552f7358SJed Brown 
1405552f7358SJed Brown   PetscFunctionBegin;
1406552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
14079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
14089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &num_proc));
14099566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
14109566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMPLEX));
1411a5b23f4aSJose E. Roman   /* Open EXODUS II file and read basic information on rank 0, then broadcast to all processors */
1412dd400576SPatrick Sanan   if (rank == 0) {
14139566063dSJacob Faibussowitsch     PetscCall(PetscMemzero(title,PETSC_MAX_PATH_LEN+1));
1414a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_init,exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs);
141528b400f6SJacob Faibussowitsch     PetscCheck(num_cs,PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set");
1416552f7358SJed Brown   }
14179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm));
14189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&dim, 1, MPI_INT, 0, comm));
14199566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dm, title));
14209566063dSJacob Faibussowitsch   PetscCall(DMPlexSetChart(*dm, 0, numCells+numVertices));
1421412e9a14SMatthew G. Knepley   /*   We do not want this label automatically computed, instead we compute it here */
14229566063dSJacob Faibussowitsch   PetscCall(DMCreateLabel(*dm, "celltype"));
1423552f7358SJed Brown 
1424552f7358SJed Brown   /* Read cell sets information */
1425dd400576SPatrick Sanan   if (rank == 0) {
1426552f7358SJed Brown     PetscInt *cone;
1427e8f6893fSMatthew G. Knepley     int      c, cs, ncs, c_loc, v, v_loc;
1428552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
1429e8f6893fSMatthew G. Knepley     int *cs_id, *cs_order;
1430552f7358SJed Brown     /* Read from ex_get_elem_block() */
1431552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
1432e8f6893fSMatthew G. Knepley     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1433552f7358SJed Brown     /* Read from ex_get_elem_conn() */
1434552f7358SJed Brown     int *cs_connect;
1435552f7358SJed Brown 
1436552f7358SJed Brown     /* Get cell sets IDs */
14379566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order));
1438a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_ids,exoid, EX_ELEM_BLOCK, cs_id);
1439552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
1440552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1441e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
1442e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
14438f861fbcSMatthew G. Knepley       DMPolytopeType ct;
14448f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
14458f861fbcSMatthew G. Knepley 
14469566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1447a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type);
14489566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
14498f861fbcSMatthew G. Knepley       dim  = PetscMax(dim, DMPolytopeTypeGetDim(ct));
1450a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr);
14518f861fbcSMatthew G. Knepley       switch (ct) {
14528f861fbcSMatthew G. Knepley         case DM_POLYTOPE_TRI_PRISM:
1453e8f6893fSMatthew G. Knepley           cs_order[cs] = cs;
1454e8f6893fSMatthew G. Knepley           ++num_hybrid;
1455e8f6893fSMatthew G. Knepley           break;
1456e8f6893fSMatthew G. Knepley         default:
1457e8f6893fSMatthew G. Knepley           for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
1458e8f6893fSMatthew G. Knepley           cs_order[cs-num_hybrid] = cs;
1459e8f6893fSMatthew G. Knepley       }
1460e8f6893fSMatthew G. Knepley     }
1461552f7358SJed Brown     /* First set sizes */
1462e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
14638f861fbcSMatthew G. Knepley       DMPolytopeType ct;
14648f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
1465e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
14668f861fbcSMatthew G. Knepley 
14679566063dSJacob Faibussowitsch       PetscCall(PetscArrayzero(elem_type, sizeof(elem_type)));
1468a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_elem_type,exoid, cs_id[cs], elem_type);
14699566063dSJacob Faibussowitsch       PetscCall(ExodusGetCellType_Internal(elem_type, &ct));
1470a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr);
1471552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
14729566063dSJacob Faibussowitsch         PetscCall(DMPlexSetConeSize(*dm, c, num_vertex_per_cell));
14739566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCellType(*dm, c, ct));
1474552f7358SJed Brown       }
1475552f7358SJed Brown     }
14769566063dSJacob Faibussowitsch     for (v = numCells; v < numCells+numVertices; ++v) PetscCall(DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT));
14779566063dSJacob Faibussowitsch     PetscCall(DMSetUp(*dm));
1478e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1479e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
1480a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_block,exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr);
14819566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone));
1482a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_conn,exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL);
1483eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1484552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1485636e64ffSStefano Zampini         DMPolytopeType ct;
1486636e64ffSStefano Zampini 
1487552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
1488552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
1489552f7358SJed Brown         }
14909566063dSJacob Faibussowitsch         PetscCall(DMPlexGetCellType(*dm, c, &ct));
14919566063dSJacob Faibussowitsch         PetscCall(DMPlexInvertCell(ct, cone));
14929566063dSJacob Faibussowitsch         PetscCall(DMPlexSetCone(*dm, c, cone));
14939566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]));
1494552f7358SJed Brown       }
14959566063dSJacob Faibussowitsch       PetscCall(PetscFree2(cs_connect,cone));
1496552f7358SJed Brown     }
14979566063dSJacob Faibussowitsch     PetscCall(PetscFree2(cs_id, cs_order));
1498552f7358SJed Brown   }
14998f861fbcSMatthew G. Knepley   {
15008f861fbcSMatthew G. Knepley     PetscInt ints[] = {dim, dimEmbed};
15018f861fbcSMatthew G. Knepley 
15029566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(ints, 2, MPIU_INT, 0, comm));
15039566063dSJacob Faibussowitsch     PetscCall(DMSetDimension(*dm, ints[0]));
15049566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinateDim(*dm, ints[1]));
15058f861fbcSMatthew G. Knepley     dim      = ints[0];
15068f861fbcSMatthew G. Knepley     dimEmbed = ints[1];
15078f861fbcSMatthew G. Knepley   }
15089566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(*dm));
15099566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(*dm));
1510552f7358SJed Brown   if (interpolate) {
15115fd9971aSMatthew G. Knepley     DM idm;
1512552f7358SJed Brown 
15139566063dSJacob Faibussowitsch     PetscCall(DMPlexInterpolate(*dm, &idm));
15149566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dm));
1515552f7358SJed Brown     *dm  = idm;
1516552f7358SJed Brown   }
1517552f7358SJed Brown 
1518552f7358SJed Brown   /* Create vertex set label */
1519dd400576SPatrick Sanan   if (rank == 0 && (num_vs > 0)) {
1520552f7358SJed Brown     int vs, v;
1521552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1522552f7358SJed Brown     int *vs_id;
1523552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1524f958083aSBlaise Bourdin     int num_vertex_in_set;
1525552f7358SJed Brown     /* Read from ex_get_node_set() */
1526552f7358SJed Brown     int *vs_vertex_list;
1527552f7358SJed Brown 
1528552f7358SJed Brown     /* Get vertex set ids */
15299566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_vs, &vs_id));
1530a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_ids,exoid, EX_NODE_SET, vs_id);
1531552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
1532a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_set_param,exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL);
15339566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(num_vertex_in_set, &vs_vertex_list));
1534a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_set,exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL);
1535552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
15369566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]));
1537552f7358SJed Brown       }
15389566063dSJacob Faibussowitsch       PetscCall(PetscFree(vs_vertex_list));
1539552f7358SJed Brown     }
15409566063dSJacob Faibussowitsch     PetscCall(PetscFree(vs_id));
1541552f7358SJed Brown   }
1542552f7358SJed Brown   /* Read coordinates */
15439566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(*dm, &coordSection));
15449566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetNumFields(coordSection, 1));
15459566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetFieldComponents(coordSection, 0, dimEmbed));
15469566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(coordSection, numCells, numCells + numVertices));
1547552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
15489566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetDof(coordSection, v, dimEmbed));
15499566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed));
1550552f7358SJed Brown   }
15519566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(coordSection));
15529566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(coordSection, &coordSize));
15539566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &coordinates));
15549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) coordinates, "coordinates"));
15559566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(coordinates, coordSize, PETSC_DETERMINE));
15569566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(coordinates, dimEmbed));
15579566063dSJacob Faibussowitsch   PetscCall(VecSetType(coordinates,VECSTANDARD));
15589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(coordinates, &coords));
1559dd400576SPatrick Sanan   if (rank == 0) {
1560e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1561552f7358SJed Brown 
15629566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z));
1563a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_coord,exoid, x, y, z);
15648f861fbcSMatthew G. Knepley     if (dimEmbed > 0) {
15658f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+0] = x[v];
15660d644c17SKarl Rupp     }
15678f861fbcSMatthew G. Knepley     if (dimEmbed > 1) {
15688f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+1] = y[v];
15690d644c17SKarl Rupp     }
15708f861fbcSMatthew G. Knepley     if (dimEmbed > 2) {
15718f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+2] = z[v];
15720d644c17SKarl Rupp     }
15739566063dSJacob Faibussowitsch     PetscCall(PetscFree3(x,y,z));
1574552f7358SJed Brown   }
15759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(coordinates, &coords));
15769566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinatesLocal(*dm, coordinates));
15779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&coordinates));
1578552f7358SJed Brown 
1579552f7358SJed Brown   /* Create side set label */
1580dd400576SPatrick Sanan   if (rank == 0 && interpolate && (num_fs > 0)) {
1581552f7358SJed Brown     int fs, f, voff;
1582552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1583552f7358SJed Brown     int *fs_id;
1584552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1585f958083aSBlaise Bourdin     int num_side_in_set;
1586552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1587552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
1588ef073283Smichael_afanasiev     /* Read side set labels */
1589ef073283Smichael_afanasiev     char fs_name[MAX_STR_LENGTH+1];
1590552f7358SJed Brown 
1591552f7358SJed Brown     /* Get side set ids */
15929566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(num_fs, &fs_id));
1593a74df02fSJacob Faibussowitsch     PetscStackCallStandard(ex_get_ids,exoid, EX_SIDE_SET, fs_id);
1594552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
1595a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_set_param,exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL);
15969566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list));
1597a74df02fSJacob Faibussowitsch       PetscStackCallStandard(ex_get_side_set_node_list,exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list);
1598ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1599ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1600552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
16010298fd71SBarry Smith         const PetscInt *faces   = NULL;
1602552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1603552f7358SJed Brown         PetscInt       faceVertices[4], v;
1604552f7358SJed Brown 
16052c71b3e2SJacob Faibussowitsch         PetscCheckFalse(faceSize > 4,comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1606552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
1607552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1608552f7358SJed Brown         }
16099566063dSJacob Faibussowitsch         PetscCall(DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
16102c71b3e2SJacob Faibussowitsch         PetscCheckFalse(numFaces != 1,comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
16119566063dSJacob Faibussowitsch         PetscCall(DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]));
1612ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1613ef073283Smichael_afanasiev         if (!fs_name_err) {
16149566063dSJacob Faibussowitsch           PetscCall(DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]));
1615ef073283Smichael_afanasiev         }
16169566063dSJacob Faibussowitsch         PetscCall(DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces));
1617552f7358SJed Brown       }
16189566063dSJacob Faibussowitsch       PetscCall(PetscFree2(fs_vertex_count_list,fs_vertex_list));
1619552f7358SJed Brown     }
16209566063dSJacob Faibussowitsch     PetscCall(PetscFree(fs_id));
1621552f7358SJed Brown   }
1622ae9eebabSLisandro Dalcin 
1623ae9eebabSLisandro Dalcin   { /* Create Cell/Face/Vertex Sets labels at all processes */
1624ae9eebabSLisandro Dalcin     enum {n = 3};
1625ae9eebabSLisandro Dalcin     PetscBool flag[n];
1626ae9eebabSLisandro Dalcin 
1627ae9eebabSLisandro Dalcin     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1628ae9eebabSLisandro Dalcin     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1629ae9eebabSLisandro Dalcin     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
16309566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(flag, n, MPIU_BOOL, 0, comm));
16319566063dSJacob Faibussowitsch     if (flag[0]) PetscCall(DMCreateLabel(*dm, "Cell Sets"));
16329566063dSJacob Faibussowitsch     if (flag[1]) PetscCall(DMCreateLabel(*dm, "Face Sets"));
16339566063dSJacob Faibussowitsch     if (flag[2]) PetscCall(DMCreateLabel(*dm, "Vertex Sets"));
1634ae9eebabSLisandro Dalcin   }
1635b458e8f1SJose E. Roman   PetscFunctionReturn(0);
1636552f7358SJed Brown #else
1637552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1638552f7358SJed Brown #endif
1639552f7358SJed Brown }
1640