xref: /petsc/src/dm/impls/plex/plexexodusii.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
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   PetscErrorCode        ierr;
491e50132fSMatthew G. Knepley 
501e50132fSMatthew G. Knepley   PetscFunctionBegin;
511e50132fSMatthew G. Knepley   if (exo->filename) {ierr = PetscViewerASCIIPrintf(viewer, "Filename:    %s\n", exo->filename);CHKERRQ(ierr);}
526823f3c5SBlaise Bourdin   if (exo->exoid)    {ierr = PetscViewerASCIIPrintf(viewer, "exoid:       %d\n", exo->exoid)   ;CHKERRQ(ierr);}
536823f3c5SBlaise Bourdin   if (exo->btype)    {ierr = PetscViewerASCIIPrintf(viewer, "IO Mode:     %d\n", exo->btype)   ;CHKERRQ(ierr);}
546823f3c5SBlaise Bourdin   if (exo->order)    {ierr = PetscViewerASCIIPrintf(viewer, "Mesh order:  %d\n", exo->order)   ;CHKERRQ(ierr);}
551e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
561e50132fSMatthew G. Knepley }
571e50132fSMatthew G. Knepley 
581e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetFromOptions_ExodusII(PetscOptionItems *PetscOptionsObject, PetscViewer v)
591e50132fSMatthew G. Knepley {
601e50132fSMatthew G. Knepley   PetscErrorCode ierr;
611e50132fSMatthew G. Knepley 
621e50132fSMatthew G. Knepley   PetscFunctionBegin;
631e50132fSMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "ExodusII PetscViewer Options");CHKERRQ(ierr);
641e50132fSMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
651e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
661e50132fSMatthew G. Knepley }
671e50132fSMatthew G. Knepley 
681e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerSetUp_ExodusII(PetscViewer viewer)
691e50132fSMatthew G. Knepley {
701e50132fSMatthew G. Knepley   PetscFunctionBegin;
711e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
721e50132fSMatthew G. Knepley }
731e50132fSMatthew G. Knepley 
741e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerDestroy_ExodusII(PetscViewer viewer)
751e50132fSMatthew G. Knepley {
761e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
771e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
781e50132fSMatthew G. Knepley 
791e50132fSMatthew G. Knepley   PetscFunctionBegin;
806823f3c5SBlaise Bourdin   if (exo->exoid >= 0) {PetscStackCallStandard(ex_close,(exo->exoid));}
816823f3c5SBlaise Bourdin   ierr = PetscFree(exo->filename);CHKERRQ(ierr);
821e50132fSMatthew G. Knepley   ierr = PetscFree(exo);CHKERRQ(ierr);
831e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
841e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
851e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
866823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetId_C",NULL);CHKERRQ(ierr);
876823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerGetOrder_C",NULL);CHKERRQ(ierr);
886823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerSetOrder_C",NULL);CHKERRQ(ierr);
891e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
901e50132fSMatthew G. Knepley }
911e50132fSMatthew G. Knepley 
921e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetName_ExodusII(PetscViewer viewer, const char name[])
931e50132fSMatthew G. Knepley {
941e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
951e50132fSMatthew G. Knepley   PetscMPIInt           rank;
961e50132fSMatthew G. Knepley   int                   CPU_word_size, IO_word_size, EXO_mode;
971e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
986823f3c5SBlaise Bourdin   MPI_Info              mpi_info = MPI_INFO_NULL;
996823f3c5SBlaise Bourdin   float                 EXO_version;
1001e50132fSMatthew G. Knepley 
101ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) viewer), &rank);CHKERRMPI(ierr);
1021e50132fSMatthew G. Knepley   CPU_word_size = sizeof(PetscReal);
1031e50132fSMatthew G. Knepley   IO_word_size  = sizeof(PetscReal);
1041e50132fSMatthew G. Knepley 
1051e50132fSMatthew G. Knepley   PetscFunctionBegin;
1066823f3c5SBlaise Bourdin   if (exo->exoid >= 0) {
1076823f3c5SBlaise Bourdin     PetscStackCallStandard(ex_close,(exo->exoid));
1086823f3c5SBlaise Bourdin     exo->exoid = -1;
1096823f3c5SBlaise Bourdin   }
1101e50132fSMatthew G. Knepley   if (exo->filename) {ierr = PetscFree(exo->filename);CHKERRQ(ierr);}
1111e50132fSMatthew G. Knepley   ierr = PetscStrallocpy(name, &exo->filename);CHKERRQ(ierr);
1121e50132fSMatthew G. Knepley   switch (exo->btype) {
1131e50132fSMatthew G. Knepley   case FILE_MODE_READ:
1146823f3c5SBlaise Bourdin     EXO_mode = EX_READ;
1151e50132fSMatthew G. Knepley     break;
1161e50132fSMatthew G. Knepley   case FILE_MODE_APPEND:
1176823f3c5SBlaise Bourdin   case FILE_MODE_UPDATE:
1186823f3c5SBlaise Bourdin   case FILE_MODE_APPEND_UPDATE:
1196823f3c5SBlaise Bourdin     /* Will fail if the file does not already exist */
1206823f3c5SBlaise Bourdin     EXO_mode = EX_WRITE;
1211e50132fSMatthew G. Knepley     break;
1221e50132fSMatthew G. Knepley   case FILE_MODE_WRITE:
1236823f3c5SBlaise Bourdin     /*
1246823f3c5SBlaise Bourdin       exodus only allows writing geometry upon file creation, so we will let DMView create the file.
1256823f3c5SBlaise Bourdin     */
1266823f3c5SBlaise Bourdin     PetscFunctionReturn(0);
1271e50132fSMatthew G. Knepley     break;
1281e50132fSMatthew G. Knepley   default:
1291e50132fSMatthew G. Knepley     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
1301e50132fSMatthew G. Knepley   }
1311e50132fSMatthew G. Knepley   #if defined(PETSC_USE_64BIT_INDICES)
1321e50132fSMatthew G. Knepley   EXO_mode += EX_ALL_INT64_API;
1331e50132fSMatthew G. Knepley   #endif
1346823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(name,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,PETSC_COMM_WORLD,mpi_info);
1356823f3c5SBlaise Bourdin   if (exo->exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", name);
1361e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1371e50132fSMatthew G. Knepley }
1381e50132fSMatthew G. Knepley 
1391e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetName_ExodusII(PetscViewer viewer, const char **name)
1401e50132fSMatthew G. Knepley {
1411e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1421e50132fSMatthew G. Knepley 
1431e50132fSMatthew G. Knepley   PetscFunctionBegin;
1441e50132fSMatthew G. Knepley   *name = exo->filename;
1451e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1461e50132fSMatthew G. Knepley }
1471e50132fSMatthew G. Knepley 
1481e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileSetMode_ExodusII(PetscViewer viewer, PetscFileMode type)
1491e50132fSMatthew G. Knepley {
1501e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1511e50132fSMatthew G. Knepley 
1521e50132fSMatthew G. Knepley   PetscFunctionBegin;
1531e50132fSMatthew G. Knepley   exo->btype = type;
1541e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1551e50132fSMatthew G. Knepley }
1561e50132fSMatthew G. Knepley 
1571e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerFileGetMode_ExodusII(PetscViewer viewer, PetscFileMode *type)
1581e50132fSMatthew G. Knepley {
1591e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1601e50132fSMatthew G. Knepley 
1611e50132fSMatthew G. Knepley   PetscFunctionBegin;
1621e50132fSMatthew G. Knepley   *type = exo->btype;
1631e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1641e50132fSMatthew G. Knepley }
1651e50132fSMatthew G. Knepley 
1661e50132fSMatthew G. Knepley static PetscErrorCode PetscViewerExodusIIGetId_ExodusII(PetscViewer viewer, int *exoid)
1671e50132fSMatthew G. Knepley {
1681e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1691e50132fSMatthew G. Knepley 
1701e50132fSMatthew G. Knepley   PetscFunctionBegin;
1711e50132fSMatthew G. Knepley   *exoid = exo->exoid;
1721e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1731e50132fSMatthew G. Knepley }
1741e50132fSMatthew G. Knepley 
1756823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIIGetOrder_ExodusII(PetscViewer viewer, PetscInt *order)
1761e50132fSMatthew G. Knepley {
1776823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1781e50132fSMatthew G. Knepley 
1791e50132fSMatthew G. Knepley   PetscFunctionBegin;
1806823f3c5SBlaise Bourdin   *order = exo->order;
1816823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
1826823f3c5SBlaise Bourdin }
1836823f3c5SBlaise Bourdin 
1846823f3c5SBlaise Bourdin static PetscErrorCode PetscViewerExodusIISetOrder_ExodusII(PetscViewer viewer, PetscInt order)
1856823f3c5SBlaise Bourdin {
1866823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
1876823f3c5SBlaise Bourdin 
1886823f3c5SBlaise Bourdin   PetscFunctionBegin;
1896823f3c5SBlaise Bourdin   exo->order = order;
1901e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
1911e50132fSMatthew G. Knepley }
1921e50132fSMatthew G. Knepley 
1931e50132fSMatthew G. Knepley /*MC
19400373969SVaclav Hapla    PETSCVIEWEREXODUSII - A viewer that writes to an Exodus II file
1951e50132fSMatthew G. Knepley 
19600373969SVaclav Hapla .seealso:  PetscViewerExodusIIOpen(), PetscViewerCreate(), PETSCVIEWERBINARY, PETSCVIEWERHDF5, DMView(),
1971e50132fSMatthew G. Knepley            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
1981e50132fSMatthew G. Knepley 
1991e50132fSMatthew G. Knepley   Level: beginner
2001e50132fSMatthew G. Knepley M*/
2011e50132fSMatthew G. Knepley 
2021e50132fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerCreate_ExodusII(PetscViewer v)
2031e50132fSMatthew G. Knepley {
2041e50132fSMatthew G. Knepley   PetscViewer_ExodusII *exo;
2051e50132fSMatthew G. Knepley   PetscErrorCode        ierr;
2061e50132fSMatthew G. Knepley 
2071e50132fSMatthew G. Knepley   PetscFunctionBegin;
2081e50132fSMatthew G. Knepley   ierr = PetscNewLog(v,&exo);CHKERRQ(ierr);
2091e50132fSMatthew G. Knepley 
2101e50132fSMatthew G. Knepley   v->data                = (void*) exo;
2111e50132fSMatthew G. Knepley   v->ops->destroy        = PetscViewerDestroy_ExodusII;
2121e50132fSMatthew G. Knepley   v->ops->setfromoptions = PetscViewerSetFromOptions_ExodusII;
2131e50132fSMatthew G. Knepley   v->ops->setup          = PetscViewerSetUp_ExodusII;
2141e50132fSMatthew G. Knepley   v->ops->view           = PetscViewerView_ExodusII;
2151e50132fSMatthew G. Knepley   v->ops->flush          = 0;
2161e50132fSMatthew G. Knepley   exo->btype             = (PetscFileMode) -1;
2171e50132fSMatthew G. Knepley   exo->filename          = 0;
2181e50132fSMatthew G. Knepley   exo->exoid             = -1;
2191e50132fSMatthew G. Knepley 
2201e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_ExodusII);CHKERRQ(ierr);
2211e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_ExodusII);CHKERRQ(ierr);
2221e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_ExodusII);CHKERRQ(ierr);
2231e50132fSMatthew G. Knepley   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_ExodusII);CHKERRQ(ierr);
2246823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetId_C",PetscViewerExodusIIGetId_ExodusII);CHKERRQ(ierr);
2256823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerSetOrder_C",PetscViewerExodusIISetOrder_ExodusII);CHKERRQ(ierr);
2266823f3c5SBlaise Bourdin   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerGetOrder_C",PetscViewerExodusIIGetOrder_ExodusII);CHKERRQ(ierr);
2271e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
2281e50132fSMatthew G. Knepley }
2291e50132fSMatthew G. Knepley 
2301e50132fSMatthew G. Knepley /*
23178569bb4SMatthew G. Knepley   EXOGetVarIndex - Locate a result in an exodus file based on its name
23278569bb4SMatthew G. Knepley 
2336823f3c5SBlaise Bourdin   Collective
23478569bb4SMatthew G. Knepley 
23578569bb4SMatthew G. Knepley   Input Parameters:
23678569bb4SMatthew G. Knepley + exoid    - the exodus id of a file (obtained from ex_open or ex_create for instance)
23778569bb4SMatthew G. Knepley . obj_type - the type of entity for instance EX_NODAL, EX_ELEM_BLOCK
23878569bb4SMatthew G. Knepley - name     - the name of the result
23978569bb4SMatthew G. Knepley 
24078569bb4SMatthew G. Knepley   Output Parameters:
24178569bb4SMatthew G. Knepley . varIndex - the location in the exodus file of the result
24278569bb4SMatthew G. Knepley 
24378569bb4SMatthew G. Knepley   Notes:
24478569bb4SMatthew G. Knepley   The exodus variable index is obtained by comparing name and the
24578569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance if name is "V"
24678569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
24778569bb4SMatthew G. Knepley   amongst all variables of type obj_type.
24878569bb4SMatthew G. Knepley 
24978569bb4SMatthew G. Knepley   Level: beginner
25078569bb4SMatthew G. Knepley 
251cd921712SStefano Zampini .seealso: EXOGetVarIndex(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadNodal_PlexEXO(),VecLoadZonal_PlexEXO()
25278569bb4SMatthew G. Knepley */
2536823f3c5SBlaise Bourdin PetscErrorCode EXOGetVarIndex_Internal(int exoid, ex_entity_type obj_type, const char name[], int *varIndex)
25478569bb4SMatthew G. Knepley {
2556de834b4SMatthew G. Knepley   int            num_vars, i, j;
25678569bb4SMatthew G. Knepley   char           ext_name[MAX_STR_LENGTH+1], var_name[MAX_STR_LENGTH+1];
25778569bb4SMatthew G. Knepley   const int      num_suffix = 5;
25878569bb4SMatthew G. Knepley   char          *suffix[5];
25978569bb4SMatthew G. Knepley   PetscBool      flg;
26078569bb4SMatthew G. Knepley   PetscErrorCode ierr;
26178569bb4SMatthew G. Knepley 
26278569bb4SMatthew G. Knepley   PetscFunctionBegin;
26378569bb4SMatthew G. Knepley   suffix[0] = (char *) "";
26478569bb4SMatthew G. Knepley   suffix[1] = (char *) "_X";
26578569bb4SMatthew G. Knepley   suffix[2] = (char *) "_XX";
26678569bb4SMatthew G. Knepley   suffix[3] = (char *) "_1";
26778569bb4SMatthew G. Knepley   suffix[4] = (char *) "_11";
26878569bb4SMatthew G. Knepley 
2696823f3c5SBlaise Bourdin   *varIndex = -1;
2706de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_variable_param,(exoid, obj_type, &num_vars));
27178569bb4SMatthew G. Knepley   for (i = 0; i < num_vars; ++i) {
2726de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_variable_name,(exoid, obj_type, i+1, var_name));
27378569bb4SMatthew G. Knepley     for (j = 0; j < num_suffix; ++j) {
27478569bb4SMatthew G. Knepley       ierr = PetscStrncpy(ext_name, name, MAX_STR_LENGTH);CHKERRQ(ierr);
275a126751eSBarry Smith       ierr = PetscStrlcat(ext_name, suffix[j], MAX_STR_LENGTH);CHKERRQ(ierr);
276*1e1ea65dSPierre Jolivet       ierr = PetscStrcasecmp(ext_name, var_name, &flg);CHKERRQ(ierr);
27778569bb4SMatthew G. Knepley       if (flg) {
27878569bb4SMatthew G. Knepley         *varIndex = i+1;
2796823f3c5SBlaise Bourdin       }
2806823f3c5SBlaise Bourdin     }
2816823f3c5SBlaise Bourdin   }
28278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
28378569bb4SMatthew G. Knepley }
28478569bb4SMatthew G. Knepley 
28578569bb4SMatthew G. Knepley /*
2866823f3c5SBlaise Bourdin   DMView_PlexExodusII - Write a DM to disk in exodus format
28778569bb4SMatthew G. Knepley 
28878569bb4SMatthew G. Knepley   Collective on dm
28978569bb4SMatthew G. Knepley 
29078569bb4SMatthew G. Knepley   Input Parameters:
29178569bb4SMatthew G. Knepley + dm  - The dm to be written
2926823f3c5SBlaise Bourdin . viewer - an exodusII viewer
29378569bb4SMatthew G. Knepley 
29478569bb4SMatthew G. Knepley   Notes:
29578569bb4SMatthew G. Knepley   Not all DM can be written to disk this way. For instance, exodus assume that element blocks (mapped to "Cell sets" labels)
29678569bb4SMatthew G. Knepley   consists of sequentially numbered cells. If this is not the case, the exodus file will be corrupted.
29778569bb4SMatthew G. Knepley 
2986823f3c5SBlaise Bourdin   If the dm has been distributed, only the part of the DM on MPI rank 0 (including "ghost" cells and vertices)
2996823f3c5SBlaise Bourdin   will be written.
30078569bb4SMatthew G. Knepley 
30178569bb4SMatthew G. Knepley   DMPlex only represents geometry while most post-processing software expect that a mesh also provides information
30278569bb4SMatthew G. Knepley   on the discretization space. This function assumes that the file represents Lagrange finite elements of order 1 or 2.
3036823f3c5SBlaise Bourdin   The order of the mesh shall be set using PetscViewerExodusIISetOrder
30478569bb4SMatthew G. Knepley   It should be extended to use PetscFE objects.
30578569bb4SMatthew G. Knepley 
3066823f3c5SBlaise Bourdin   This function will only handle TRI, TET, QUAD, and HEX cells.
30778569bb4SMatthew G. Knepley   Level: beginner
30878569bb4SMatthew G. Knepley 
3096823f3c5SBlaise Bourdin .seealso:
31078569bb4SMatthew G. Knepley */
3116823f3c5SBlaise Bourdin PetscErrorCode DMView_PlexExodusII(DM dm, PetscViewer viewer)
31278569bb4SMatthew G. Knepley {
31378569bb4SMatthew G. Knepley   enum ElemType {TRI, QUAD, TET, HEX};
31478569bb4SMatthew G. Knepley   MPI_Comm        comm;
3156823f3c5SBlaise Bourdin   PetscInt        degree; /* the order of the mesh */
31678569bb4SMatthew G. Knepley   /* Connectivity Variables */
31778569bb4SMatthew G. Knepley   PetscInt        cellsNotInConnectivity;
31878569bb4SMatthew G. Knepley   /* Cell Sets */
31978569bb4SMatthew G. Knepley   DMLabel         csLabel;
32078569bb4SMatthew G. Knepley   IS              csIS;
32178569bb4SMatthew G. Knepley   const PetscInt *csIdx;
32278569bb4SMatthew G. Knepley   PetscInt        num_cs, cs;
32378569bb4SMatthew G. Knepley   enum ElemType  *type;
324fe209ef9SBlaise Bourdin   PetscBool       hasLabel;
32578569bb4SMatthew G. Knepley   /* Coordinate Variables */
32678569bb4SMatthew G. Knepley   DM              cdm;
3276823f3c5SBlaise Bourdin   PetscSection    coordSection;
32878569bb4SMatthew G. Knepley   Vec             coord;
32978569bb4SMatthew G. Knepley   PetscInt      **nodes;
330eae8a387SMatthew G. Knepley   PetscInt        depth, d, dim, skipCells = 0;
33178569bb4SMatthew G. Knepley   PetscInt        pStart, pEnd, p, cStart, cEnd, numCells, vStart, vEnd, numVertices, eStart, eEnd, numEdges, fStart, fEnd, numFaces, numNodes;
33278569bb4SMatthew G. Knepley   PetscInt        num_vs, num_fs;
33378569bb4SMatthew G. Knepley   PetscMPIInt     rank, size;
33478569bb4SMatthew G. Knepley   const char     *dmName;
335fe209ef9SBlaise Bourdin   PetscInt        nodesTriP1[4]  = {3,0,0,0};
336fe209ef9SBlaise Bourdin   PetscInt        nodesTriP2[4]  = {3,3,0,0};
337fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP1[4] = {4,0,0,0};
338fe209ef9SBlaise Bourdin   PetscInt        nodesQuadP2[4] = {4,4,0,1};
339fe209ef9SBlaise Bourdin   PetscInt        nodesTetP1[4]  = {4,0,0,0};
340fe209ef9SBlaise Bourdin   PetscInt        nodesTetP2[4]  = {4,6,0,0};
341fe209ef9SBlaise Bourdin   PetscInt        nodesHexP1[4]  = {8,0,0,0};
342fe209ef9SBlaise Bourdin   PetscInt        nodesHexP2[4]  = {8,12,6,1};
34378569bb4SMatthew G. Knepley   PetscErrorCode  ierr;
34478569bb4SMatthew G. Knepley 
3456823f3c5SBlaise Bourdin   int             CPU_word_size, IO_word_size, EXO_mode;
3466823f3c5SBlaise Bourdin   float           EXO_version;
3476823f3c5SBlaise Bourdin 
3486823f3c5SBlaise Bourdin   PetscViewer_ExodusII *exo = (PetscViewer_ExodusII *) viewer->data;
3496823f3c5SBlaise Bourdin 
35078569bb4SMatthew G. Knepley   PetscFunctionBegin;
35178569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
352ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
353ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
3546823f3c5SBlaise Bourdin 
3556823f3c5SBlaise Bourdin   /*
3566823f3c5SBlaise Bourdin     Creating coordSection is a collective operation so we do it somewhat out of sequence
3576823f3c5SBlaise Bourdin   */
3586823f3c5SBlaise Bourdin   ierr = PetscSectionCreate(comm, &coordSection);CHKERRQ(ierr);
3596823f3c5SBlaise Bourdin   ierr = DMGetCoordinatesLocalSetUp(dm);CHKERRQ(ierr);
3606823f3c5SBlaise Bourdin   if (!rank) {
3616823f3c5SBlaise Bourdin     switch (exo->btype) {
3626823f3c5SBlaise Bourdin     case FILE_MODE_READ:
3636823f3c5SBlaise Bourdin     case FILE_MODE_APPEND:
3646823f3c5SBlaise Bourdin     case FILE_MODE_UPDATE:
3656823f3c5SBlaise Bourdin     case FILE_MODE_APPEND_UPDATE:
3666823f3c5SBlaise Bourdin       /* exodusII does not allow writing geometry to an existing file */
3676823f3c5SBlaise Bourdin       SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "cannot add geometry to existing file %s", exo->filename);
3686823f3c5SBlaise Bourdin     case FILE_MODE_WRITE:
3696823f3c5SBlaise Bourdin       /* Create an empty file if one already exists*/
3706823f3c5SBlaise Bourdin       EXO_mode = EX_CLOBBER;
3716823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
3726823f3c5SBlaise Bourdin       EXO_mode += EX_ALL_INT64_API;
3736823f3c5SBlaise Bourdin #endif
3746823f3c5SBlaise Bourdin         CPU_word_size = sizeof(PetscReal);
3756823f3c5SBlaise Bourdin         IO_word_size  = sizeof(PetscReal);
3766823f3c5SBlaise Bourdin         exo->exoid = ex_create(exo->filename, EXO_mode, &CPU_word_size, &IO_word_size);
3776823f3c5SBlaise Bourdin         if (exo->exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_create failed for %s", exo->filename);
3786823f3c5SBlaise Bourdin 
3796823f3c5SBlaise Bourdin       break;
3806823f3c5SBlaise Bourdin     default:
3816823f3c5SBlaise Bourdin       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3826823f3c5SBlaise Bourdin     }
3836823f3c5SBlaise Bourdin 
38478569bb4SMatthew G. Knepley     /* --- Get DM info --- */
38578569bb4SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) dm, &dmName);CHKERRQ(ierr);
38678569bb4SMatthew G. Knepley     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
38778569bb4SMatthew G. Knepley     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
38878569bb4SMatthew G. Knepley     ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
38978569bb4SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
39078569bb4SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr);
39178569bb4SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);
39278569bb4SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
39378569bb4SMatthew G. Knepley     numCells    = cEnd - cStart;
39478569bb4SMatthew G. Knepley     numEdges    = eEnd - eStart;
39578569bb4SMatthew G. Knepley     numVertices = vEnd - vStart;
39678569bb4SMatthew G. Knepley     if (depth == 3) {numFaces = fEnd - fStart;}
39778569bb4SMatthew G. Knepley     else            {numFaces = 0;}
39878569bb4SMatthew G. Knepley     ierr = DMGetLabelSize(dm, "Cell Sets", &num_cs);CHKERRQ(ierr);
39978569bb4SMatthew G. Knepley     ierr = DMGetLabelSize(dm, "Vertex Sets", &num_vs);CHKERRQ(ierr);
40078569bb4SMatthew G. Knepley     ierr = DMGetLabelSize(dm, "Face Sets", &num_fs);CHKERRQ(ierr);
40178569bb4SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coord);CHKERRQ(ierr);
40278569bb4SMatthew G. Knepley     ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
40378569bb4SMatthew G. Knepley     if (num_cs > 0) {
40478569bb4SMatthew G. Knepley       ierr = DMGetLabel(dm, "Cell Sets", &csLabel);CHKERRQ(ierr);
40578569bb4SMatthew G. Knepley       ierr = DMLabelGetValueIS(csLabel, &csIS);CHKERRQ(ierr);
40678569bb4SMatthew G. Knepley       ierr = ISGetIndices(csIS, &csIdx);CHKERRQ(ierr);
40778569bb4SMatthew G. Knepley     }
40878569bb4SMatthew G. Knepley     ierr = PetscMalloc1(num_cs, &nodes);CHKERRQ(ierr);
40978569bb4SMatthew G. Knepley     /* Set element type for each block and compute total number of nodes */
41078569bb4SMatthew G. Knepley     ierr = PetscMalloc1(num_cs, &type);CHKERRQ(ierr);
41178569bb4SMatthew G. Knepley     numNodes = numVertices;
4126823f3c5SBlaise Bourdin 
4136823f3c5SBlaise Bourdin     ierr = PetscViewerExodusIIGetOrder(viewer, &degree);CHKERRQ(ierr);
41478569bb4SMatthew G. Knepley     if (degree == 2) {numNodes += numEdges;}
41578569bb4SMatthew G. Knepley     cellsNotInConnectivity = numCells;
41678569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
41778569bb4SMatthew G. Knepley       IS              stratumIS;
41878569bb4SMatthew G. Knepley       const PetscInt *cells;
41978569bb4SMatthew G. Knepley       PetscScalar    *xyz = NULL;
42078569bb4SMatthew G. Knepley       PetscInt        csSize, closureSize;
42178569bb4SMatthew G. Knepley 
42278569bb4SMatthew G. Knepley       ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
42378569bb4SMatthew G. Knepley       ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
42478569bb4SMatthew G. Knepley       ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
42578569bb4SMatthew G. Knepley       ierr = DMPlexVecGetClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
42678569bb4SMatthew G. Knepley       switch (dim) {
42778569bb4SMatthew G. Knepley         case 2:
42878569bb4SMatthew G. Knepley           if      (closureSize == 3*dim) {type[cs] = TRI;}
42978569bb4SMatthew G. Knepley           else if (closureSize == 4*dim) {type[cs] = QUAD;}
43078569bb4SMatthew G. Knepley           else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
43178569bb4SMatthew G. Knepley           break;
43278569bb4SMatthew G. Knepley         case 3:
43378569bb4SMatthew G. Knepley           if      (closureSize == 4*dim) {type[cs] = TET;}
43478569bb4SMatthew G. Knepley           else if (closureSize == 8*dim) {type[cs] = HEX;}
43578569bb4SMatthew G. Knepley           else SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of vertices %D in dimension %D has no ExodusII type", closureSize/dim, dim);
43678569bb4SMatthew G. Knepley           break;
43778569bb4SMatthew G. Knepley         default: SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %D not handled by ExodusII viewer", dim);
43878569bb4SMatthew G. Knepley       }
43978569bb4SMatthew G. Knepley       if ((degree == 2) && (type[cs] == QUAD)) {numNodes += csSize;}
44078569bb4SMatthew G. Knepley       if ((degree == 2) && (type[cs] == HEX))  {numNodes += csSize; numNodes += numFaces;}
44178569bb4SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(cdm, NULL, coord, cells[0], &closureSize, &xyz);CHKERRQ(ierr);
44278569bb4SMatthew G. Knepley       /* Set nodes and Element type */
44378569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
44478569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesTriP1;
44578569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTriP2;
44678569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
44778569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesQuadP1;
44878569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesQuadP2;
44978569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
45078569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesTetP1;
45178569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesTetP2;
45278569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
45378569bb4SMatthew G. Knepley         if      (degree == 1) nodes[cs] = nodesHexP1;
45478569bb4SMatthew G. Knepley         else if (degree == 2) nodes[cs] = nodesHexP2;
45578569bb4SMatthew G. Knepley       }
45678569bb4SMatthew G. Knepley       /* Compute the number of cells not in the connectivity table */
45778569bb4SMatthew G. Knepley       cellsNotInConnectivity -= nodes[cs][3]*csSize;
45878569bb4SMatthew G. Knepley 
45978569bb4SMatthew G. Knepley       ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
46078569bb4SMatthew G. Knepley       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
46178569bb4SMatthew G. Knepley     }
4626823f3c5SBlaise Bourdin     if (num_cs > 0) {PetscStackCallStandard(ex_put_init,(exo->exoid, dmName, dim, numNodes, numCells, num_cs, num_vs, num_fs));}
46378569bb4SMatthew G. Knepley     /* --- Connectivity --- */
46478569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
46578569bb4SMatthew G. Knepley       IS              stratumIS;
46678569bb4SMatthew G. Knepley       const PetscInt *cells;
467fe209ef9SBlaise Bourdin       PetscInt       *connect, off = 0;
468eae8a387SMatthew G. Knepley       PetscInt        edgesInClosure = 0, facesInClosure = 0, verticesInClosure = 0;
46978569bb4SMatthew G. Knepley       PetscInt        csSize, c, connectSize, closureSize;
47078569bb4SMatthew G. Knepley       char           *elem_type = NULL;
47178569bb4SMatthew G. Knepley       char            elem_type_tri3[]  = "TRI3",  elem_type_quad4[] = "QUAD4";
47278569bb4SMatthew G. Knepley       char            elem_type_tri6[]  = "TRI6",  elem_type_quad9[] = "QUAD9";
47378569bb4SMatthew G. Knepley       char            elem_type_tet4[]  = "TET4",  elem_type_hex8[]  = "HEX8";
47478569bb4SMatthew G. Knepley       char            elem_type_tet10[] = "TET10", elem_type_hex27[] = "HEX27";
47578569bb4SMatthew G. Knepley 
47678569bb4SMatthew G. Knepley       ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
47778569bb4SMatthew G. Knepley       ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
47878569bb4SMatthew G. Knepley       ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
47978569bb4SMatthew G. Knepley       /* Set Element type */
48078569bb4SMatthew G. Knepley       if (type[cs] == TRI) {
48178569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_tri3;
48278569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tri6;
48378569bb4SMatthew G. Knepley       } else if (type[cs] == QUAD) {
48478569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_quad4;
48578569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_quad9;
48678569bb4SMatthew G. Knepley       } else if (type[cs] == TET) {
48778569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_tet4;
48878569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_tet10;
48978569bb4SMatthew G. Knepley       } else if (type[cs] == HEX) {
49078569bb4SMatthew G. Knepley         if      (degree == 1) elem_type = elem_type_hex8;
49178569bb4SMatthew G. Knepley         else if (degree == 2) elem_type = elem_type_hex27;
49278569bb4SMatthew G. Knepley       }
49378569bb4SMatthew G. Knepley       connectSize = nodes[cs][0] + nodes[cs][1] + nodes[cs][2] + nodes[cs][3];
494fe209ef9SBlaise Bourdin       ierr = PetscMalloc1(PetscMax(27,connectSize)*csSize, &connect);CHKERRQ(ierr);
4956823f3c5SBlaise Bourdin       PetscStackCallStandard(ex_put_block,(exo->exoid, EX_ELEM_BLOCK, csIdx[cs], elem_type, csSize, connectSize, 0, 0, 1));
49678569bb4SMatthew G. Knepley       /* Find number of vertices, edges, and faces in the closure */
49778569bb4SMatthew G. Knepley       verticesInClosure = nodes[cs][0];
49878569bb4SMatthew G. Knepley       if (depth > 1) {
49978569bb4SMatthew G. Knepley         if (dim == 2) {
50078569bb4SMatthew G. Knepley           ierr = DMPlexGetConeSize(dm, cells[0], &edgesInClosure);CHKERRQ(ierr);
50178569bb4SMatthew G. Knepley         } else if (dim == 3) {
50278569bb4SMatthew G. Knepley           PetscInt *closure = NULL;
50378569bb4SMatthew G. Knepley 
50478569bb4SMatthew G. Knepley           ierr = DMPlexGetConeSize(dm, cells[0], &facesInClosure);CHKERRQ(ierr);
50578569bb4SMatthew G. Knepley           ierr = DMPlexGetTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
50678569bb4SMatthew G. Knepley           edgesInClosure = closureSize - facesInClosure - 1 - verticesInClosure;
50778569bb4SMatthew G. Knepley           ierr = DMPlexRestoreTransitiveClosure(dm, cells[0], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
50878569bb4SMatthew G. Knepley         }
50978569bb4SMatthew G. Knepley       }
51078569bb4SMatthew G. Knepley       /* Get connectivity for each cell */
51178569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
51278569bb4SMatthew G. Knepley         PetscInt *closure = NULL;
51378569bb4SMatthew G. Knepley         PetscInt  temp, i;
51478569bb4SMatthew G. Knepley 
51578569bb4SMatthew G. Knepley         ierr = DMPlexGetTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
51678569bb4SMatthew G. Knepley         for (i = 0; i < connectSize; ++i) {
51778569bb4SMatthew G. Knepley           if (i < nodes[cs][0]) {/* Vertices */
518fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i+edgesInClosure+facesInClosure+1)*2] + 1;
519fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
52078569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]) { /* Edges */
521fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i-verticesInClosure+facesInClosure+1)*2] + 1;
522fe209ef9SBlaise Bourdin             if (nodes[cs][2] == 0) connect[i+off] -= numFaces;
523fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
52478569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]) { /* Cells */
525fe209ef9SBlaise Bourdin             connect[i+off] = closure[0] + 1;
526fe209ef9SBlaise Bourdin             connect[i+off] -= skipCells;
52778569bb4SMatthew G. Knepley           } else if (i < nodes[cs][0]+nodes[cs][1]+nodes[cs][3]+nodes[cs][2]) { /* Faces */
528fe209ef9SBlaise Bourdin             connect[i+off] = closure[(i-edgesInClosure-verticesInClosure)*2] + 1;
529fe209ef9SBlaise Bourdin             connect[i+off] -= cellsNotInConnectivity;
53078569bb4SMatthew G. Knepley           } else {
531fe209ef9SBlaise Bourdin             connect[i+off] = -1;
53278569bb4SMatthew G. Knepley           }
53378569bb4SMatthew G. Knepley         }
53478569bb4SMatthew G. Knepley         /* Tetrahedra are inverted */
53578569bb4SMatthew G. Knepley         if (type[cs] == TET) {
536fe209ef9SBlaise Bourdin           temp = connect[0+off]; connect[0+off] = connect[1+off]; connect[1+off] = temp;
53778569bb4SMatthew G. Knepley           if (degree == 2) {
538e2c9586dSBlaise Bourdin             temp = connect[5+off]; connect[5+off] = connect[6+off]; connect[6+off] = temp;
539e2c9586dSBlaise Bourdin             temp = connect[7+off]; connect[7+off] = connect[8+off]; connect[8+off] = temp;
54078569bb4SMatthew G. Knepley           }
54178569bb4SMatthew G. Knepley         }
54278569bb4SMatthew G. Knepley         /* Hexahedra are inverted */
54378569bb4SMatthew G. Knepley         if (type[cs] == HEX) {
544fe209ef9SBlaise Bourdin           temp = connect[1+off]; connect[1+off] = connect[3+off]; connect[3+off] = temp;
54578569bb4SMatthew G. Knepley           if (degree == 2) {
546fe209ef9SBlaise Bourdin             temp = connect[8+off];  connect[8+off]  = connect[11+off]; connect[11+off] = temp;
547fe209ef9SBlaise Bourdin             temp = connect[9+off];  connect[9+off]  = connect[10+off]; connect[10+off] = temp;
548fe209ef9SBlaise Bourdin             temp = connect[16+off]; connect[16+off] = connect[17+off]; connect[17+off] = temp;
549fe209ef9SBlaise Bourdin             temp = connect[18+off]; connect[18+off] = connect[19+off]; connect[19+off] = temp;
55078569bb4SMatthew G. Knepley 
551fe209ef9SBlaise Bourdin             temp = connect[12+off]; connect[12+off] = connect[16+off]; connect[16+off] = temp;
552fe209ef9SBlaise Bourdin             temp = connect[13+off]; connect[13+off] = connect[17+off]; connect[17+off] = temp;
553fe209ef9SBlaise Bourdin             temp = connect[14+off]; connect[14+off] = connect[18+off]; connect[18+off] = temp;
554fe209ef9SBlaise Bourdin             temp = connect[15+off]; connect[15+off] = connect[19+off]; connect[19+off] = temp;
55578569bb4SMatthew G. Knepley 
556fe209ef9SBlaise Bourdin             temp = connect[23+off]; connect[23+off] = connect[26+off]; connect[26+off] = temp;
557fe209ef9SBlaise Bourdin             temp = connect[24+off]; connect[24+off] = connect[25+off]; connect[25+off] = temp;
558fe209ef9SBlaise Bourdin             temp = connect[25+off]; connect[25+off] = connect[26+off]; connect[26+off] = temp;
55978569bb4SMatthew G. Knepley           }
56078569bb4SMatthew G. Knepley         }
561fe209ef9SBlaise Bourdin         off += connectSize;
56278569bb4SMatthew G. Knepley         ierr = DMPlexRestoreTransitiveClosure(dm, cells[c], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
56378569bb4SMatthew G. Knepley       }
5646823f3c5SBlaise Bourdin       PetscStackCallStandard(ex_put_conn,(exo->exoid, EX_ELEM_BLOCK, csIdx[cs], connect, 0, 0));
56578569bb4SMatthew G. Knepley       skipCells += (nodes[cs][3] == 0)*csSize;
56678569bb4SMatthew G. Knepley       ierr = PetscFree(connect);CHKERRQ(ierr);
56778569bb4SMatthew G. Knepley       ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
56878569bb4SMatthew G. Knepley       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
56978569bb4SMatthew G. Knepley     }
57078569bb4SMatthew G. Knepley     ierr = PetscFree(type);CHKERRQ(ierr);
57178569bb4SMatthew G. Knepley     /* --- Coordinates --- */
5726823f3c5SBlaise Bourdin     ierr = PetscSectionSetChart(coordSection, pStart, pEnd);CHKERRQ(ierr);
5731e50132fSMatthew G. Knepley     if (num_cs) {
57478569bb4SMatthew G. Knepley       for (d = 0; d < depth; ++d) {
57578569bb4SMatthew G. Knepley         ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
57678569bb4SMatthew G. Knepley         for (p = pStart; p < pEnd; ++p) {
5776823f3c5SBlaise Bourdin           ierr = PetscSectionSetDof(coordSection, p, nodes[0][d] > 0);CHKERRQ(ierr);
57878569bb4SMatthew G. Knepley         }
57978569bb4SMatthew G. Knepley       }
5801e50132fSMatthew G. Knepley     }
58178569bb4SMatthew G. Knepley     for (cs = 0; cs < num_cs; ++cs) {
58278569bb4SMatthew G. Knepley       IS              stratumIS;
58378569bb4SMatthew G. Knepley       const PetscInt *cells;
58478569bb4SMatthew G. Knepley       PetscInt        csSize, c;
58578569bb4SMatthew G. Knepley 
58678569bb4SMatthew G. Knepley       ierr = DMLabelGetStratumIS(csLabel, csIdx[cs], &stratumIS);CHKERRQ(ierr);
58778569bb4SMatthew G. Knepley       ierr = ISGetIndices(stratumIS, &cells);CHKERRQ(ierr);
58878569bb4SMatthew G. Knepley       ierr = ISGetSize(stratumIS, &csSize);CHKERRQ(ierr);
58978569bb4SMatthew G. Knepley       for (c = 0; c < csSize; ++c) {
5906823f3c5SBlaise Bourdin         ierr = PetscSectionSetDof(coordSection, cells[c], nodes[cs][3] > 0);CHKERRQ(ierr);
59178569bb4SMatthew G. Knepley       }
59278569bb4SMatthew G. Knepley       ierr = ISRestoreIndices(stratumIS, &cells);CHKERRQ(ierr);
59378569bb4SMatthew G. Knepley       ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
59478569bb4SMatthew G. Knepley     }
59578569bb4SMatthew G. Knepley     if (num_cs > 0) {
59678569bb4SMatthew G. Knepley       ierr = ISRestoreIndices(csIS, &csIdx);CHKERRQ(ierr);
59778569bb4SMatthew G. Knepley       ierr = ISDestroy(&csIS);CHKERRQ(ierr);
59878569bb4SMatthew G. Knepley     }
59978569bb4SMatthew G. Knepley     ierr = PetscFree(nodes);CHKERRQ(ierr);
6006823f3c5SBlaise Bourdin     ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
60178569bb4SMatthew G. Knepley     if (numNodes > 0) {
60278569bb4SMatthew G. Knepley       const char  *coordNames[3] = {"x", "y", "z"};
60378569bb4SMatthew G. Knepley       PetscScalar *coords, *closure;
60478569bb4SMatthew G. Knepley       PetscReal   *cval;
60578569bb4SMatthew G. Knepley       PetscInt     hasDof, n = 0;
60678569bb4SMatthew G. Knepley 
6076823f3c5SBlaise Bourdin       /* There can't be more than 24 values in the closure of a point for the coord coordSection */
6086de834b4SMatthew G. Knepley       ierr = PetscCalloc3(numNodes*3, &coords, dim, &cval, 24, &closure);CHKERRQ(ierr);
6096823f3c5SBlaise Bourdin       ierr = DMGetCoordinatesLocalNoncollective(dm, &coord);CHKERRQ(ierr);
61078569bb4SMatthew G. Knepley       ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
61178569bb4SMatthew G. Knepley       for (p = pStart; p < pEnd; ++p) {
612*1e1ea65dSPierre Jolivet         ierr = PetscSectionGetDof(coordSection, p, &hasDof);CHKERRQ(ierr);
61378569bb4SMatthew G. Knepley         if (hasDof) {
61478569bb4SMatthew G. Knepley           PetscInt closureSize = 24, j;
61578569bb4SMatthew G. Knepley 
61678569bb4SMatthew G. Knepley           ierr = DMPlexVecGetClosure(cdm, NULL, coord, p, &closureSize, &closure);CHKERRQ(ierr);
61778569bb4SMatthew G. Knepley           for (d = 0; d < dim; ++d) {
61878569bb4SMatthew G. Knepley             cval[d] = 0.0;
61978569bb4SMatthew G. Knepley             for (j = 0; j < closureSize/dim; j++) cval[d] += closure[j*dim+d];
62078569bb4SMatthew G. Knepley             coords[d*numNodes+n] = cval[d] * dim / closureSize;
62178569bb4SMatthew G. Knepley           }
62278569bb4SMatthew G. Knepley           ++n;
62378569bb4SMatthew G. Knepley         }
62478569bb4SMatthew G. Knepley       }
6256823f3c5SBlaise Bourdin       PetscStackCallStandard(ex_put_coord,(exo->exoid, &coords[0*numNodes], &coords[1*numNodes], &coords[2*numNodes]));
62678569bb4SMatthew G. Knepley       ierr = PetscFree3(coords, cval, closure);CHKERRQ(ierr);
6276823f3c5SBlaise Bourdin       PetscStackCallStandard(ex_put_coord_names,(exo->exoid, (char **) coordNames));
62878569bb4SMatthew G. Knepley     }
6296823f3c5SBlaise Bourdin 
630fe209ef9SBlaise Bourdin     /* --- Node Sets/Vertex Sets --- */
631fe209ef9SBlaise Bourdin     DMHasLabel(dm, "Vertex Sets", &hasLabel);
632fe209ef9SBlaise Bourdin     if (hasLabel) {
633fe209ef9SBlaise Bourdin       PetscInt        i, vs, vsSize;
634fe209ef9SBlaise Bourdin       const PetscInt *vsIdx, *vertices;
635fe209ef9SBlaise Bourdin       PetscInt       *nodeList;
636fe209ef9SBlaise Bourdin       IS              vsIS, stratumIS;
637fe209ef9SBlaise Bourdin       DMLabel         vsLabel;
638fe209ef9SBlaise Bourdin       ierr = DMGetLabel(dm, "Vertex Sets", &vsLabel);CHKERRQ(ierr);
639fe209ef9SBlaise Bourdin       ierr = DMLabelGetValueIS(vsLabel, &vsIS);CHKERRQ(ierr);
640fe209ef9SBlaise Bourdin       ierr = ISGetIndices(vsIS, &vsIdx);CHKERRQ(ierr);
641fe209ef9SBlaise Bourdin       for (vs=0; vs<num_vs; ++vs) {
642fe209ef9SBlaise Bourdin         ierr = DMLabelGetStratumIS(vsLabel, vsIdx[vs], &stratumIS);CHKERRQ(ierr);
643fe209ef9SBlaise Bourdin         ierr = ISGetIndices(stratumIS, &vertices);CHKERRQ(ierr);
644fe209ef9SBlaise Bourdin         ierr = ISGetSize(stratumIS, &vsSize);CHKERRQ(ierr);
645*1e1ea65dSPierre Jolivet         ierr = PetscMalloc1(vsSize, &nodeList);CHKERRQ(ierr);
646fe209ef9SBlaise Bourdin         for (i=0; i<vsSize; ++i) {
647fe209ef9SBlaise Bourdin           nodeList[i] = vertices[i] - skipCells + 1;
648fe209ef9SBlaise Bourdin         }
6496823f3c5SBlaise Bourdin         PetscStackCallStandard(ex_put_set_param,(exo->exoid, EX_NODE_SET, vsIdx[vs], vsSize, 0));
6506823f3c5SBlaise Bourdin         PetscStackCallStandard(ex_put_set,(exo->exoid, EX_NODE_SET, vsIdx[vs], nodeList, NULL));
651fe209ef9SBlaise Bourdin         ierr = ISRestoreIndices(stratumIS, &vertices);CHKERRQ(ierr);
652fe209ef9SBlaise Bourdin         ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
653fe209ef9SBlaise Bourdin         ierr = PetscFree(nodeList);CHKERRQ(ierr);
654fe209ef9SBlaise Bourdin       }
655fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(vsIS, &vsIdx);CHKERRQ(ierr);
656fe209ef9SBlaise Bourdin       ierr = ISDestroy(&vsIS);CHKERRQ(ierr);
657fe209ef9SBlaise Bourdin     }
658fe209ef9SBlaise Bourdin     /* --- Side Sets/Face Sets --- */
659fe209ef9SBlaise Bourdin     ierr = DMHasLabel(dm, "Face Sets", &hasLabel);CHKERRQ(ierr);
660fe209ef9SBlaise Bourdin     if (hasLabel) {
661fe209ef9SBlaise Bourdin       PetscInt        i, j, fs, fsSize;
662fe209ef9SBlaise Bourdin       const PetscInt *fsIdx, *faces;
663fe209ef9SBlaise Bourdin       IS              fsIS, stratumIS;
664fe209ef9SBlaise Bourdin       DMLabel         fsLabel;
665fe209ef9SBlaise Bourdin       PetscInt        numPoints, *points;
666fe209ef9SBlaise Bourdin       PetscInt        elem_list_size = 0;
667fe209ef9SBlaise Bourdin       PetscInt       *elem_list, *elem_ind, *side_list;
668fe209ef9SBlaise Bourdin 
669fe209ef9SBlaise Bourdin       ierr = DMGetLabel(dm, "Face Sets", &fsLabel);CHKERRQ(ierr);
670fe209ef9SBlaise Bourdin       /* Compute size of Node List and Element List */
671fe209ef9SBlaise Bourdin       ierr = DMLabelGetValueIS(fsLabel, &fsIS);CHKERRQ(ierr);
672fe209ef9SBlaise Bourdin       ierr = ISGetIndices(fsIS, &fsIdx);CHKERRQ(ierr);
673fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
674feb237baSPierre Jolivet         ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);
675fe209ef9SBlaise Bourdin         ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
676fe209ef9SBlaise Bourdin         elem_list_size += fsSize;
677fe209ef9SBlaise Bourdin         ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
678fe209ef9SBlaise Bourdin       }
679fe209ef9SBlaise Bourdin       ierr = PetscMalloc3(num_fs, &elem_ind,
680fe209ef9SBlaise Bourdin                           elem_list_size, &elem_list,
681fe209ef9SBlaise Bourdin                           elem_list_size, &side_list);CHKERRQ(ierr);
682fe209ef9SBlaise Bourdin       elem_ind[0] = 0;
683fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
684fe209ef9SBlaise Bourdin         ierr = DMLabelGetStratumIS(fsLabel, fsIdx[fs], &stratumIS);CHKERRQ(ierr);
685fe209ef9SBlaise Bourdin         ierr = ISGetIndices(stratumIS, &faces);CHKERRQ(ierr);
686fe209ef9SBlaise Bourdin         ierr = ISGetSize(stratumIS, &fsSize);CHKERRQ(ierr);
687fe209ef9SBlaise Bourdin         /* Set Parameters */
6886823f3c5SBlaise Bourdin         PetscStackCallStandard(ex_put_set_param,(exo->exoid, EX_SIDE_SET, fsIdx[fs], fsSize, 0));
689fe209ef9SBlaise Bourdin         /* Indices */
690fe209ef9SBlaise Bourdin         if (fs<num_fs-1) {
691fe209ef9SBlaise Bourdin           elem_ind[fs+1] = elem_ind[fs] + fsSize;
692fe209ef9SBlaise Bourdin         }
693fe209ef9SBlaise Bourdin 
694fe209ef9SBlaise Bourdin         for (i=0; i<fsSize; ++i) {
695fe209ef9SBlaise Bourdin           /* Element List */
696fe209ef9SBlaise Bourdin           points = NULL;
697fe209ef9SBlaise Bourdin           ierr = DMPlexGetTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
698fe209ef9SBlaise Bourdin           elem_list[elem_ind[fs] + i] = points[2] +1;
699fe209ef9SBlaise Bourdin           ierr = DMPlexRestoreTransitiveClosure(dm, faces[i], PETSC_FALSE, &numPoints, &points);CHKERRQ(ierr);
700fe209ef9SBlaise Bourdin 
701fe209ef9SBlaise Bourdin           /* Side List */
702fe209ef9SBlaise Bourdin           points = NULL;
703fe209ef9SBlaise Bourdin           ierr = DMPlexGetTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
704fe209ef9SBlaise Bourdin           for (j=1; j<numPoints; ++j) {
705fe209ef9SBlaise Bourdin             if (points[j*2]==faces[i]) {break;}
706fe209ef9SBlaise Bourdin           }
707fe209ef9SBlaise Bourdin           /* Convert HEX sides */
708fe209ef9SBlaise Bourdin           if (numPoints == 27) {
709fe209ef9SBlaise Bourdin             if      (j == 1) {j = 5;}
710fe209ef9SBlaise Bourdin             else if (j == 2) {j = 6;}
711fe209ef9SBlaise Bourdin             else if (j == 3) {j = 1;}
712fe209ef9SBlaise Bourdin             else if (j == 4) {j = 3;}
713fe209ef9SBlaise Bourdin             else if (j == 5) {j = 2;}
714fe209ef9SBlaise Bourdin             else if (j == 6) {j = 4;}
715fe209ef9SBlaise Bourdin           }
716fe209ef9SBlaise Bourdin           /* Convert TET sides */
717fe209ef9SBlaise Bourdin           if (numPoints == 15) {
718fe209ef9SBlaise Bourdin             --j;
719fe209ef9SBlaise Bourdin             if (j == 0) {j = 4;}
720fe209ef9SBlaise Bourdin           }
721fe209ef9SBlaise Bourdin           side_list[elem_ind[fs] + i] = j;
722fe209ef9SBlaise Bourdin           ierr = DMPlexRestoreTransitiveClosure(dm, elem_list[elem_ind[fs] + i]-1, PETSC_TRUE, &numPoints, &points);CHKERRQ(ierr);
723fe209ef9SBlaise Bourdin 
724fe209ef9SBlaise Bourdin         }
725fe209ef9SBlaise Bourdin         ierr = ISRestoreIndices(stratumIS, &faces);CHKERRQ(ierr);
726fe209ef9SBlaise Bourdin         ierr = ISDestroy(&stratumIS);CHKERRQ(ierr);
727fe209ef9SBlaise Bourdin       }
728fe209ef9SBlaise Bourdin       ierr = ISRestoreIndices(fsIS, &fsIdx);CHKERRQ(ierr);
729fe209ef9SBlaise Bourdin       ierr = ISDestroy(&fsIS);CHKERRQ(ierr);
730fe209ef9SBlaise Bourdin 
731fe209ef9SBlaise Bourdin       /* Put side sets */
732fe209ef9SBlaise Bourdin       for (fs=0; fs<num_fs; ++fs) {
7336823f3c5SBlaise Bourdin         PetscStackCallStandard(ex_put_set,(exo->exoid, EX_SIDE_SET, fsIdx[fs], &elem_list[elem_ind[fs]], &side_list[elem_ind[fs]]));
734fe209ef9SBlaise Bourdin       }
735fe209ef9SBlaise Bourdin       ierr = PetscFree3(elem_ind,elem_list,side_list);CHKERRQ(ierr);
736fe209ef9SBlaise Bourdin     }
7376823f3c5SBlaise Bourdin     /*
7386823f3c5SBlaise Bourdin       close the exodus file
7396823f3c5SBlaise Bourdin     */
7406823f3c5SBlaise Bourdin     ex_close(exo->exoid);
7416823f3c5SBlaise Bourdin     exo->exoid = -1;
7426823f3c5SBlaise Bourdin   }
7436823f3c5SBlaise Bourdin   ierr = PetscSectionDestroy(&coordSection);CHKERRQ(ierr);
7446823f3c5SBlaise Bourdin 
7456823f3c5SBlaise Bourdin   /*
7466823f3c5SBlaise Bourdin     reopen the file in parallel
7476823f3c5SBlaise Bourdin   */
7486823f3c5SBlaise Bourdin   EXO_mode = EX_WRITE;
7496823f3c5SBlaise Bourdin #if defined(PETSC_USE_64BIT_INDICES)
7506823f3c5SBlaise Bourdin   EXO_mode += EX_ALL_INT64_API;
7516823f3c5SBlaise Bourdin #endif
7526823f3c5SBlaise Bourdin   CPU_word_size = sizeof(PetscReal);
7536823f3c5SBlaise Bourdin   IO_word_size  = sizeof(PetscReal);
7546823f3c5SBlaise Bourdin   exo->exoid = ex_open_par(exo->filename,EXO_mode,&CPU_word_size,&IO_word_size,&EXO_version,comm,MPI_INFO_NULL);
7556823f3c5SBlaise Bourdin   if (exo->exoid < 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open_par failed for %s", exo->filename);
7566823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
7576823f3c5SBlaise Bourdin }
7586823f3c5SBlaise Bourdin 
7596823f3c5SBlaise Bourdin /*
7606823f3c5SBlaise Bourdin   VecView_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
7616823f3c5SBlaise Bourdin 
7626823f3c5SBlaise Bourdin   Collective on v
7636823f3c5SBlaise Bourdin 
7646823f3c5SBlaise Bourdin   Input Parameters:
7656823f3c5SBlaise Bourdin + v  - The vector to be written
7666823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
7676823f3c5SBlaise Bourdin 
7686823f3c5SBlaise Bourdin   Notes:
7696823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
7706823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
7716823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
7726823f3c5SBlaise Bourdin   amongst all variables.
7736823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
7746823f3c5SBlaise Bourdin   possibly corrupting the file
7756823f3c5SBlaise Bourdin 
7766823f3c5SBlaise Bourdin   Level: beginner
7776823f3c5SBlaise Bourdin 
7786823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII(),VecView_PlexExodusII()
7796823f3c5SBlaise Bourdin @*/
7806823f3c5SBlaise Bourdin PetscErrorCode VecView_PlexExodusII_Internal(Vec v, PetscViewer viewer)
7816823f3c5SBlaise Bourdin {
7826823f3c5SBlaise Bourdin   DM                 dm;
7836823f3c5SBlaise Bourdin   MPI_Comm           comm;
7846823f3c5SBlaise Bourdin   PetscMPIInt        rank;
7856823f3c5SBlaise Bourdin 
7866823f3c5SBlaise Bourdin   int                exoid,offsetN = 0, offsetZ = 0;
7876823f3c5SBlaise Bourdin   const char        *vecname;
7886823f3c5SBlaise Bourdin   PetscInt           step;
7896823f3c5SBlaise Bourdin   PetscErrorCode     ierr;
7906823f3c5SBlaise Bourdin 
7916823f3c5SBlaise Bourdin   PetscFunctionBegin;
7926823f3c5SBlaise Bourdin   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
7936823f3c5SBlaise Bourdin   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
7946823f3c5SBlaise Bourdin   ierr = PetscViewerExodusIIGetId(viewer,&exoid);CHKERRQ(ierr);
7956823f3c5SBlaise Bourdin   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
7966823f3c5SBlaise Bourdin   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
7976823f3c5SBlaise Bourdin 
7986823f3c5SBlaise Bourdin   ierr = DMGetOutputSequenceNumber(dm,&step,NULL);CHKERRQ(ierr);
7996823f3c5SBlaise Bourdin   ierr = EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN);CHKERRQ(ierr);
8006823f3c5SBlaise Bourdin   ierr = EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ);CHKERRQ(ierr);
8016823f3c5SBlaise Bourdin   if (offsetN <= 0 && offsetZ <= 0) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. \n", vecname);
8026823f3c5SBlaise Bourdin   if (offsetN > 0) {
8036823f3c5SBlaise Bourdin     ierr = VecViewPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN);CHKERRQ(ierr);
8046823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
8056823f3c5SBlaise Bourdin     ierr = VecViewPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ);CHKERRQ(ierr);
8066823f3c5SBlaise Bourdin   } else SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. \n", vecname);
8076823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
8086823f3c5SBlaise Bourdin }
8096823f3c5SBlaise Bourdin 
8106823f3c5SBlaise Bourdin /*
8116823f3c5SBlaise Bourdin   VecLoad_PlexExodusII_Internal - Write a Vec corresponding in an exodus file
8126823f3c5SBlaise Bourdin 
8136823f3c5SBlaise Bourdin   Collective on v
8146823f3c5SBlaise Bourdin 
8156823f3c5SBlaise Bourdin   Input Parameters:
8166823f3c5SBlaise Bourdin + v  - The vector to be written
8176823f3c5SBlaise Bourdin - viewer - The PetscViewerExodusII viewer associate to an exodus file
8186823f3c5SBlaise Bourdin 
8196823f3c5SBlaise Bourdin   Notes:
8206823f3c5SBlaise Bourdin   The exodus result variable index is obtained by comparing the Vec name and the
8216823f3c5SBlaise Bourdin   names of variables declared in the exodus file. For instance for a Vec named "V"
8226823f3c5SBlaise Bourdin   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
8236823f3c5SBlaise Bourdin   amongst all variables.
8246823f3c5SBlaise Bourdin   In the event where a nodal and zonal variable both match, the function will return an error instead of
8256823f3c5SBlaise Bourdin   possibly corrupting the file
8266823f3c5SBlaise Bourdin 
8276823f3c5SBlaise Bourdin   Level: beginner
8286823f3c5SBlaise Bourdin 
8296823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII(),VecView_PlexExodusII()
8306823f3c5SBlaise Bourdin @*/
8316823f3c5SBlaise Bourdin PetscErrorCode VecLoad_PlexExodusII_Internal(Vec v, PetscViewer viewer)
8326823f3c5SBlaise Bourdin {
8336823f3c5SBlaise Bourdin   DM                 dm;
8346823f3c5SBlaise Bourdin   MPI_Comm           comm;
8356823f3c5SBlaise Bourdin   PetscMPIInt        rank;
8366823f3c5SBlaise Bourdin 
8376823f3c5SBlaise Bourdin   int                exoid,offsetN = 0, offsetZ = 0;
8386823f3c5SBlaise Bourdin   const char        *vecname;
8396823f3c5SBlaise Bourdin   PetscInt           step;
8406823f3c5SBlaise Bourdin   PetscErrorCode     ierr;
8416823f3c5SBlaise Bourdin 
8426823f3c5SBlaise Bourdin   PetscFunctionBegin;
8436823f3c5SBlaise Bourdin   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
8446823f3c5SBlaise Bourdin   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
8456823f3c5SBlaise Bourdin   ierr = PetscViewerExodusIIGetId(viewer,&exoid);CHKERRQ(ierr);
8466823f3c5SBlaise Bourdin   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
8476823f3c5SBlaise Bourdin   ierr = PetscObjectGetName((PetscObject) v, &vecname);CHKERRQ(ierr);
8486823f3c5SBlaise Bourdin 
8496823f3c5SBlaise Bourdin   ierr = DMGetOutputSequenceNumber(dm,&step,NULL);CHKERRQ(ierr);
8506823f3c5SBlaise Bourdin   ierr = EXOGetVarIndex_Internal(exoid,EX_NODAL,vecname,&offsetN);CHKERRQ(ierr);
8516823f3c5SBlaise Bourdin   ierr = EXOGetVarIndex_Internal(exoid,EX_ELEM_BLOCK,vecname,&offsetZ);CHKERRQ(ierr);
8526823f3c5SBlaise Bourdin   if (offsetN <= 0 && offsetZ <= 0) SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Found both nodal and zonal variable %s in exodus file. \n", vecname);
8536823f3c5SBlaise Bourdin   if (offsetN > 0) {
8546823f3c5SBlaise Bourdin     ierr = VecLoadPlex_ExodusII_Nodal_Internal(v,exoid,(int) step+1,offsetN);CHKERRQ(ierr);
8556823f3c5SBlaise Bourdin   } else if (offsetZ > 0) {
8566823f3c5SBlaise Bourdin     ierr = VecLoadPlex_ExodusII_Zonal_Internal(v,exoid,(int) step+1,offsetZ);CHKERRQ(ierr);
8576823f3c5SBlaise Bourdin   } else SETERRQ1(comm, PETSC_ERR_FILE_UNEXPECTED, "Could not find nodal or zonal variable %s in exodus file. \n", vecname);
85878569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
85978569bb4SMatthew G. Knepley }
86078569bb4SMatthew G. Knepley 
86178569bb4SMatthew G. Knepley /*
86278569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Nodal_Internal - Write a Vec corresponding to a nodal field to an exodus file
86378569bb4SMatthew G. Knepley 
86478569bb4SMatthew G. Knepley   Collective on v
86578569bb4SMatthew G. Knepley 
86678569bb4SMatthew G. Knepley   Input Parameters:
86778569bb4SMatthew G. Knepley + v  - The vector to be written
86878569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
8696823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
8706823f3c5SBlaise Bourdin - offset - the location of the variable in the file
87178569bb4SMatthew G. Knepley 
87278569bb4SMatthew G. Knepley   Notes:
87378569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
87478569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
87578569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
87678569bb4SMatthew G. Knepley   amongst all nodal variables.
87778569bb4SMatthew G. Knepley 
87878569bb4SMatthew G. Knepley   Level: beginner
87978569bb4SMatthew G. Knepley 
8806823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII_Internal(),VecLoadNodal_PlexEXO(),VecViewZonal_PlexEXO(),VecLoadZonal_PlexEXO()
88178569bb4SMatthew G. Knepley @*/
8826823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
88378569bb4SMatthew G. Knepley {
88478569bb4SMatthew G. Knepley   MPI_Comm           comm;
88578569bb4SMatthew G. Knepley   PetscMPIInt        size;
88678569bb4SMatthew G. Knepley   DM                 dm;
88778569bb4SMatthew G. Knepley   Vec                vNatural, vComp;
88822a7544eSVaclav Hapla   const PetscScalar *varray;
88978569bb4SMatthew G. Knepley   PetscInt           xs, xe, bs;
89078569bb4SMatthew G. Knepley   PetscBool          useNatural;
89178569bb4SMatthew G. Knepley   PetscErrorCode     ierr;
89278569bb4SMatthew G. Knepley 
89378569bb4SMatthew G. Knepley   PetscFunctionBegin;
89478569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
895ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
89678569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
89778569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
89878569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
89978569bb4SMatthew G. Knepley   if (useNatural) {
90078569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
90178569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
90278569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
90378569bb4SMatthew G. Knepley   } else {
90478569bb4SMatthew G. Knepley     vNatural = v;
90578569bb4SMatthew G. Knepley   }
9066823f3c5SBlaise Bourdin 
90778569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
90878569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
90978569bb4SMatthew G. Knepley      We assume that they are stored sequentially */
91078569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
91178569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
91278569bb4SMatthew G. Knepley   if (bs == 1) {
91378569bb4SMatthew G. Knepley     ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
9146de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
91578569bb4SMatthew G. Knepley     ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
91678569bb4SMatthew G. Knepley   } else {
91778569bb4SMatthew G. Knepley     IS       compIS;
91878569bb4SMatthew G. Knepley     PetscInt c;
91978569bb4SMatthew G. Knepley 
92094a2d30dSStefano Zampini     ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
92178569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
92278569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
92378569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
92478569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
9256de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
92678569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
92778569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
92878569bb4SMatthew G. Knepley     }
92978569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
93078569bb4SMatthew G. Knepley   }
93178569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);}
93278569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
93378569bb4SMatthew G. Knepley }
93478569bb4SMatthew G. Knepley 
93578569bb4SMatthew G. Knepley /*
93678569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Nodal_Internal - Read a Vec corresponding to a nodal field from an exodus file
93778569bb4SMatthew G. Knepley 
93878569bb4SMatthew G. Knepley   Collective on v
93978569bb4SMatthew G. Knepley 
94078569bb4SMatthew G. Knepley   Input Parameters:
94178569bb4SMatthew G. Knepley + v  - The vector to be written
94278569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
9436823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
9446823f3c5SBlaise Bourdin - offset - the location of the variable in the file
94578569bb4SMatthew G. Knepley 
94678569bb4SMatthew G. Knepley   Notes:
94778569bb4SMatthew G. Knepley   The exodus result nodal variable index is obtained by comparing the Vec name and the
94878569bb4SMatthew G. Knepley   names of nodal variables declared in the exodus file. For instance for a Vec named "V"
94978569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
95078569bb4SMatthew G. Knepley   amongst all nodal variables.
95178569bb4SMatthew G. Knepley 
95278569bb4SMatthew G. Knepley   Level: beginner
95378569bb4SMatthew G. Knepley 
9546823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecViewZonal_PlexEXO(), VecLoadZonal_PlexEXO()
95578569bb4SMatthew G. Knepley */
9566823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Nodal_Internal(Vec v, int exoid, int step, int offset)
95778569bb4SMatthew G. Knepley {
95878569bb4SMatthew G. Knepley   MPI_Comm       comm;
95978569bb4SMatthew G. Knepley   PetscMPIInt    size;
96078569bb4SMatthew G. Knepley   DM             dm;
96178569bb4SMatthew G. Knepley   Vec            vNatural, vComp;
96222a7544eSVaclav Hapla   PetscScalar   *varray;
96378569bb4SMatthew G. Knepley   PetscInt       xs, xe, bs;
96478569bb4SMatthew G. Knepley   PetscBool      useNatural;
96578569bb4SMatthew G. Knepley   PetscErrorCode ierr;
96678569bb4SMatthew G. Knepley 
96778569bb4SMatthew G. Knepley   PetscFunctionBegin;
96878569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr);
969ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
97078569bb4SMatthew G. Knepley   ierr = VecGetDM(v,&dm);CHKERRQ(ierr);
97178569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
97278569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
97378569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
97478569bb4SMatthew G. Knepley   else            {vNatural = v;}
9756823f3c5SBlaise Bourdin 
97678569bb4SMatthew G. Knepley   /* Read local chunk from the file */
97778569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
97878569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
97978569bb4SMatthew G. Knepley   if (bs == 1) {
98078569bb4SMatthew G. Knepley     ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
9816de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset, 1, xs+1, xe-xs, varray));
98278569bb4SMatthew G. Knepley     ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
98378569bb4SMatthew G. Knepley   } else {
98478569bb4SMatthew G. Knepley     IS       compIS;
98578569bb4SMatthew G. Knepley     PetscInt c;
98678569bb4SMatthew G. Knepley 
98794a2d30dSStefano Zampini     ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);
98878569bb4SMatthew G. Knepley     for (c = 0; c < bs; ++c) {
98978569bb4SMatthew G. Knepley       ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
99078569bb4SMatthew G. Knepley       ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
99178569bb4SMatthew G. Knepley       ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
9926de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_NODAL, offset+c, 1, xs/bs+1, (xe-xs)/bs, varray));
99378569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
99478569bb4SMatthew G. Knepley       ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
99578569bb4SMatthew G. Knepley     }
99678569bb4SMatthew G. Knepley     ierr = ISDestroy(&compIS);CHKERRQ(ierr);
99778569bb4SMatthew G. Knepley   }
99878569bb4SMatthew G. Knepley   if (useNatural) {
99978569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
100078569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
100178569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
100278569bb4SMatthew G. Knepley   }
100378569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
100478569bb4SMatthew G. Knepley }
100578569bb4SMatthew G. Knepley 
100678569bb4SMatthew G. Knepley /*
100778569bb4SMatthew G. Knepley   VecViewPlex_ExodusII_Zonal_Internal - Write a Vec corresponding to a zonal (cell based) field to an exodus file
100878569bb4SMatthew G. Knepley 
100978569bb4SMatthew G. Knepley   Collective on v
101078569bb4SMatthew G. Knepley 
101178569bb4SMatthew G. Knepley   Input Parameters:
101278569bb4SMatthew G. Knepley + v  - The vector to be written
101378569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
10146823f3c5SBlaise Bourdin . step - the time step to write at (exodus steps are numbered starting from 1)
10156823f3c5SBlaise Bourdin - offset - the location of the variable in the file
101678569bb4SMatthew G. Knepley 
101778569bb4SMatthew G. Knepley   Notes:
101878569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
101978569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
102078569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
102178569bb4SMatthew G. Knepley   amongst all zonal variables.
102278569bb4SMatthew G. Knepley 
102378569bb4SMatthew G. Knepley   Level: beginner
102478569bb4SMatthew G. Knepley 
10256823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(),DMPlexView_ExodusII_Internal(),VecViewPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Nodal_Internal(),VecLoadPlex_ExodusII_Zonal_Internal()
102678569bb4SMatthew G. Knepley */
10276823f3c5SBlaise Bourdin PetscErrorCode VecViewPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
102878569bb4SMatthew G. Knepley {
102978569bb4SMatthew G. Knepley   MPI_Comm          comm;
103078569bb4SMatthew G. Knepley   PetscMPIInt       size;
103178569bb4SMatthew G. Knepley   DM                dm;
103278569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
103322a7544eSVaclav Hapla   const PetscScalar *varray;
103478569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
103578569bb4SMatthew G. Knepley   PetscBool         useNatural;
103678569bb4SMatthew G. Knepley   IS                compIS;
103778569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
103878569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
103978569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
104078569bb4SMatthew G. Knepley 
104178569bb4SMatthew G. Knepley   PetscFunctionBegin;
104278569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v, &comm);CHKERRQ(ierr);
1043ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
104478569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
104578569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
104678569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
104778569bb4SMatthew G. Knepley   if (useNatural) {
104878569bb4SMatthew G. Knepley     ierr = DMGetGlobalVector(dm, &vNatural);CHKERRQ(ierr);
104978569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalBegin(dm, v, vNatural);CHKERRQ(ierr);
105078569bb4SMatthew G. Knepley     ierr = DMPlexGlobalToNaturalEnd(dm, v, vNatural);CHKERRQ(ierr);
105178569bb4SMatthew G. Knepley   } else {
105278569bb4SMatthew G. Knepley     vNatural = v;
105378569bb4SMatthew G. Knepley   }
10546823f3c5SBlaise Bourdin 
105578569bb4SMatthew G. Knepley   /* Write local chunk of the result in the exodus file
105678569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
105778569bb4SMatthew G. Knepley      We assume that they are stored sequentially
105878569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
105978569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
106078569bb4SMatthew G. Knepley      to figure out what to save where. */
106178569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
106278569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
10636de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
106478569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
106578569bb4SMatthew G. Knepley     ex_block block;
106678569bb4SMatthew G. Knepley 
106778569bb4SMatthew G. Knepley     block.id   = csID[set];
106878569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
10696de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
107078569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
107178569bb4SMatthew G. Knepley   }
107278569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
107378569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
107478569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
107578569bb4SMatthew G. Knepley   for (set = 0; set < numCS; set++) {
107678569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
107778569bb4SMatthew G. Knepley 
107878569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
107978569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
108078569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
108178569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
108278569bb4SMatthew G. Knepley     if (bs == 1) {
108378569bb4SMatthew G. Knepley       ierr = VecGetArrayRead(vNatural, &varray);CHKERRQ(ierr);
10846de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
108578569bb4SMatthew G. Knepley       ierr = VecRestoreArrayRead(vNatural, &varray);CHKERRQ(ierr);
108678569bb4SMatthew G. Knepley     } else {
108778569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
108878569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
108978569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
109078569bb4SMatthew G. Knepley         ierr = VecGetArrayRead(vComp, &varray);CHKERRQ(ierr);
10916de834b4SMatthew G. Knepley         PetscStackCallStandard(ex_put_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]));
109278569bb4SMatthew G. Knepley         ierr = VecRestoreArrayRead(vComp, &varray);CHKERRQ(ierr);
109378569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
109478569bb4SMatthew G. Knepley       }
109578569bb4SMatthew G. Knepley     }
109678569bb4SMatthew G. Knepley     csxs += csSize[set];
109778569bb4SMatthew G. Knepley   }
109878569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
109978569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
110078569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMRestoreGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
110178569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
110278569bb4SMatthew G. Knepley }
110378569bb4SMatthew G. Knepley 
110478569bb4SMatthew G. Knepley /*
110578569bb4SMatthew G. Knepley   VecLoadPlex_ExodusII_Zonal_Internal - Read a Vec corresponding to a zonal (cell based) field from an exodus file
110678569bb4SMatthew G. Knepley 
110778569bb4SMatthew G. Knepley   Collective on v
110878569bb4SMatthew G. Knepley 
110978569bb4SMatthew G. Knepley   Input Parameters:
111078569bb4SMatthew G. Knepley + v  - The vector to be written
111178569bb4SMatthew G. Knepley . exoid - the exodus id of a file (obtained from ex_open or ex_create for instance)
11126823f3c5SBlaise Bourdin . step - the time step to read at (exodus steps are numbered starting from 1)
11136823f3c5SBlaise Bourdin - offset - the location of the variable in the file
111478569bb4SMatthew G. Knepley 
111578569bb4SMatthew G. Knepley   Notes:
111678569bb4SMatthew G. Knepley   The exodus result zonal variable index is obtained by comparing the Vec name and the
111778569bb4SMatthew G. Knepley   names of zonal variables declared in the exodus file. For instance for a Vec named "V"
111878569bb4SMatthew G. Knepley   the location in the exodus file will be the first match of "V", "V_X", "V_XX", "V_1", or "V_11"
111978569bb4SMatthew G. Knepley   amongst all zonal variables.
112078569bb4SMatthew G. Knepley 
112178569bb4SMatthew G. Knepley   Level: beginner
112278569bb4SMatthew G. Knepley 
11236823f3c5SBlaise Bourdin .seealso: EXOGetVarIndex_Internal(), DMPlexView_ExodusII_Internal(), VecViewPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Nodal_Internal(), VecLoadPlex_ExodusII_Zonal_Internal()
112478569bb4SMatthew G. Knepley */
11256823f3c5SBlaise Bourdin PetscErrorCode VecLoadPlex_ExodusII_Zonal_Internal(Vec v, int exoid, int step, int offset)
112678569bb4SMatthew G. Knepley {
112778569bb4SMatthew G. Knepley   MPI_Comm          comm;
112878569bb4SMatthew G. Knepley   PetscMPIInt       size;
112978569bb4SMatthew G. Knepley   DM                dm;
113078569bb4SMatthew G. Knepley   Vec               vNatural, vComp;
113122a7544eSVaclav Hapla   PetscScalar      *varray;
113278569bb4SMatthew G. Knepley   PetscInt          xs, xe, bs;
113378569bb4SMatthew G. Knepley   PetscBool         useNatural;
113478569bb4SMatthew G. Knepley   IS                compIS;
113578569bb4SMatthew G. Knepley   PetscInt         *csSize, *csID;
113678569bb4SMatthew G. Knepley   PetscInt          numCS, set, csxs = 0;
113778569bb4SMatthew G. Knepley   PetscErrorCode    ierr;
113878569bb4SMatthew G. Knepley 
113978569bb4SMatthew G. Knepley   PetscFunctionBegin;
114078569bb4SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)v,&comm);CHKERRQ(ierr);
1141ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
114278569bb4SMatthew G. Knepley   ierr = VecGetDM(v, &dm);CHKERRQ(ierr);
114378569bb4SMatthew G. Knepley   ierr = DMGetUseNatural(dm, &useNatural);CHKERRQ(ierr);
114478569bb4SMatthew G. Knepley   useNatural = useNatural && size > 1 ? PETSC_TRUE : PETSC_FALSE;
114578569bb4SMatthew G. Knepley   if (useNatural) {ierr = DMGetGlobalVector(dm,&vNatural);CHKERRQ(ierr);}
114678569bb4SMatthew G. Knepley   else            {vNatural = v;}
11476823f3c5SBlaise Bourdin 
114878569bb4SMatthew G. Knepley   /* Read local chunk of the result in the exodus file
114978569bb4SMatthew G. Knepley      exodus stores each component of a vector-valued field as a separate variable.
115078569bb4SMatthew G. Knepley      We assume that they are stored sequentially
115178569bb4SMatthew G. Knepley      Zonal variables are accessed one element block at a time, so we loop through the cell sets,
115278569bb4SMatthew G. Knepley      but once the vector has been reordered to natural size, we cannot use the label informations
115378569bb4SMatthew G. Knepley      to figure out what to save where. */
115478569bb4SMatthew G. Knepley   numCS = ex_inquire_int(exoid, EX_INQ_ELEM_BLK);
115578569bb4SMatthew G. Knepley   ierr = PetscMalloc2(numCS, &csID, numCS, &csSize);CHKERRQ(ierr);
11566de834b4SMatthew G. Knepley   PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, csID));
115778569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
115878569bb4SMatthew G. Knepley     ex_block block;
115978569bb4SMatthew G. Knepley 
116078569bb4SMatthew G. Knepley     block.id   = csID[set];
116178569bb4SMatthew G. Knepley     block.type = EX_ELEM_BLOCK;
11626de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_block_param,(exoid, &block));
116378569bb4SMatthew G. Knepley     csSize[set] = block.num_entry;
116478569bb4SMatthew G. Knepley   }
116578569bb4SMatthew G. Knepley   ierr = VecGetOwnershipRange(vNatural, &xs, &xe);CHKERRQ(ierr);
116678569bb4SMatthew G. Knepley   ierr = VecGetBlockSize(vNatural, &bs);CHKERRQ(ierr);
116778569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISCreateStride(comm, (xe-xs)/bs, xs, bs, &compIS);CHKERRQ(ierr);}
116878569bb4SMatthew G. Knepley   for (set = 0; set < numCS; ++set) {
116978569bb4SMatthew G. Knepley     PetscInt csLocalSize, c;
117078569bb4SMatthew G. Knepley 
117178569bb4SMatthew G. Knepley     /* range of indices for set setID[set]: csxs:csxs + csSize[set]-1
117278569bb4SMatthew G. Knepley        local slice of zonal values:         xs/bs,xm/bs-1
117378569bb4SMatthew G. Knepley        intersection:                        max(xs/bs,csxs),min(xm/bs-1,csxs + csSize[set]-1) */
117478569bb4SMatthew G. Knepley     csLocalSize = PetscMax(0, PetscMin(xe/bs, csxs+csSize[set]) - PetscMax(xs/bs, csxs));
117578569bb4SMatthew G. Knepley     if (bs == 1) {
117678569bb4SMatthew G. Knepley       ierr = VecGetArray(vNatural, &varray);CHKERRQ(ierr);
11776de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset, csID[set], PetscMax(xs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs)]));
117878569bb4SMatthew G. Knepley       ierr = VecRestoreArray(vNatural, &varray);CHKERRQ(ierr);
117978569bb4SMatthew G. Knepley     } else {
118078569bb4SMatthew G. Knepley       for (c = 0; c < bs; ++c) {
118178569bb4SMatthew G. Knepley         ierr = ISStrideSetStride(compIS, (xe-xs)/bs, xs+c, bs);CHKERRQ(ierr);
118278569bb4SMatthew G. Knepley         ierr = VecGetSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
118378569bb4SMatthew G. Knepley         ierr = VecGetArray(vComp, &varray);CHKERRQ(ierr);
11846de834b4SMatthew G. Knepley         PetscStackCallStandard(ex_get_partial_var,(exoid, step, EX_ELEM_BLOCK, offset+c, csID[set], PetscMax(xs/bs-csxs, 0)+1, csLocalSize, &varray[PetscMax(0, csxs-xs/bs)]));
118578569bb4SMatthew G. Knepley         ierr = VecRestoreArray(vComp, &varray);CHKERRQ(ierr);
118678569bb4SMatthew G. Knepley         ierr = VecRestoreSubVector(vNatural, compIS, &vComp);CHKERRQ(ierr);
118778569bb4SMatthew G. Knepley       }
118878569bb4SMatthew G. Knepley     }
118978569bb4SMatthew G. Knepley     csxs += csSize[set];
119078569bb4SMatthew G. Knepley   }
119178569bb4SMatthew G. Knepley   ierr = PetscFree2(csID, csSize);CHKERRQ(ierr);
119278569bb4SMatthew G. Knepley   if (bs > 1) {ierr = ISDestroy(&compIS);CHKERRQ(ierr);}
119378569bb4SMatthew G. Knepley   if (useNatural) {
119478569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalBegin(dm, vNatural, v);CHKERRQ(ierr);
119578569bb4SMatthew G. Knepley     ierr = DMPlexNaturalToGlobalEnd(dm, vNatural, v);CHKERRQ(ierr);
119678569bb4SMatthew G. Knepley     ierr = DMRestoreGlobalVector(dm, &vNatural);CHKERRQ(ierr);
119778569bb4SMatthew G. Knepley   }
119878569bb4SMatthew G. Knepley   PetscFunctionReturn(0);
119978569bb4SMatthew G. Knepley }
1200b53c8456SSatish Balay #endif
120178569bb4SMatthew G. Knepley 
12021e50132fSMatthew G. Knepley /*@
12031e50132fSMatthew G. Knepley   PetscViewerExodusIIGetId - Get the file id of the ExodusII file
12041e50132fSMatthew G. Knepley 
12051e50132fSMatthew G. Knepley   Logically Collective on PetscViewer
12061e50132fSMatthew G. Knepley 
12071e50132fSMatthew G. Knepley   Input Parameter:
12081e50132fSMatthew G. Knepley .  viewer - the PetscViewer
12091e50132fSMatthew G. Knepley 
12101e50132fSMatthew G. Knepley   Output Parameter:
12111e50132fSMatthew G. Knepley -  exoid - The ExodusII file id
12121e50132fSMatthew G. Knepley 
12131e50132fSMatthew G. Knepley   Level: intermediate
12141e50132fSMatthew G. Knepley 
12151e50132fSMatthew G. Knepley .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
12161e50132fSMatthew G. Knepley @*/
12171e50132fSMatthew G. Knepley PetscErrorCode PetscViewerExodusIIGetId(PetscViewer viewer, int *exoid)
12181e50132fSMatthew G. Knepley {
12191e50132fSMatthew G. Knepley   PetscErrorCode ierr;
12201e50132fSMatthew G. Knepley 
12211e50132fSMatthew G. Knepley   PetscFunctionBegin;
12221e50132fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12236823f3c5SBlaise Bourdin   ierr = PetscTryMethod(viewer, "PetscViewerGetId_C",(PetscViewer,int*),(viewer,exoid));CHKERRQ(ierr);
12246823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12256823f3c5SBlaise Bourdin }
12266823f3c5SBlaise Bourdin 
12276823f3c5SBlaise Bourdin /*@
12286823f3c5SBlaise Bourdin    PetscViewerExodusIISetOrder - Set the elements order in the exodusII file.
12296823f3c5SBlaise Bourdin 
12306823f3c5SBlaise Bourdin    Collective
12316823f3c5SBlaise Bourdin 
12326823f3c5SBlaise Bourdin    Input Parameters:
12336823f3c5SBlaise Bourdin +  viewer - the viewer
123498a6a78aSPierre Jolivet -  order - elements order
12356823f3c5SBlaise Bourdin 
12366823f3c5SBlaise Bourdin    Output Parameter:
12376823f3c5SBlaise Bourdin 
12386823f3c5SBlaise Bourdin    Level: beginner
12396823f3c5SBlaise Bourdin 
12406823f3c5SBlaise Bourdin    Note:
12416823f3c5SBlaise Bourdin 
12426823f3c5SBlaise Bourdin .seealso: PetscViewerExodusIIGetId(), PetscViewerExodusIIGetOrder(), PetscViewerExodusIISetOrder()
12436823f3c5SBlaise Bourdin @*/
12446823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIISetOrder(PetscViewer viewer, PetscInt order)
12456823f3c5SBlaise Bourdin {
12466823f3c5SBlaise Bourdin   PetscErrorCode ierr;
12476823f3c5SBlaise Bourdin 
12486823f3c5SBlaise Bourdin   PetscFunctionBegin;
12496823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12506823f3c5SBlaise Bourdin   ierr = PetscTryMethod(viewer, "PetscViewerSetOrder_C",(PetscViewer,PetscInt),(viewer,order));CHKERRQ(ierr);
12516823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12526823f3c5SBlaise Bourdin }
12536823f3c5SBlaise Bourdin 
12546823f3c5SBlaise Bourdin /*@
12556823f3c5SBlaise Bourdin    PetscViewerExodusIIGetOrder - Get the elements order in the exodusII file.
12566823f3c5SBlaise Bourdin 
12576823f3c5SBlaise Bourdin    Collective
12586823f3c5SBlaise Bourdin 
12596823f3c5SBlaise Bourdin    Input Parameters:
12606823f3c5SBlaise Bourdin +  viewer - the viewer
126198a6a78aSPierre Jolivet -  order - elements order
12626823f3c5SBlaise Bourdin 
12636823f3c5SBlaise Bourdin    Output Parameter:
12646823f3c5SBlaise Bourdin 
12656823f3c5SBlaise Bourdin    Level: beginner
12666823f3c5SBlaise Bourdin 
12676823f3c5SBlaise Bourdin    Note:
12686823f3c5SBlaise Bourdin 
12696823f3c5SBlaise Bourdin .seealso: PetscViewerExodusIIGetId(), PetscViewerExodusIIGetOrder(), PetscViewerExodusIISetOrder()
12706823f3c5SBlaise Bourdin @*/
12716823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIGetOrder(PetscViewer viewer, PetscInt *order)
12726823f3c5SBlaise Bourdin {
12736823f3c5SBlaise Bourdin   PetscErrorCode ierr;
12746823f3c5SBlaise Bourdin 
12756823f3c5SBlaise Bourdin   PetscFunctionBegin;
12766823f3c5SBlaise Bourdin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12776823f3c5SBlaise Bourdin   ierr = PetscTryMethod(viewer, "PetscViewerGetOrder_C",(PetscViewer,PetscInt*),(viewer,order));CHKERRQ(ierr);
12786823f3c5SBlaise Bourdin   PetscFunctionReturn(0);
12796823f3c5SBlaise Bourdin }
12806823f3c5SBlaise Bourdin 
12816823f3c5SBlaise Bourdin /*@C
12826823f3c5SBlaise Bourdin    PetscViewerExodusIIOpen - Opens a file for ExodusII input/output.
12836823f3c5SBlaise Bourdin 
12846823f3c5SBlaise Bourdin    Collective
12856823f3c5SBlaise Bourdin 
12866823f3c5SBlaise Bourdin    Input Parameters:
12876823f3c5SBlaise Bourdin +  comm - MPI communicator
12886823f3c5SBlaise Bourdin .  name - name of file
12896823f3c5SBlaise Bourdin -  type - type of file
12906823f3c5SBlaise Bourdin $    FILE_MODE_WRITE - create new file for binary output
12916823f3c5SBlaise Bourdin $    FILE_MODE_READ - open existing file for binary input
12926823f3c5SBlaise Bourdin $    FILE_MODE_APPEND - open existing file for binary output
12936823f3c5SBlaise Bourdin 
12946823f3c5SBlaise Bourdin    Output Parameter:
12956823f3c5SBlaise Bourdin .  exo - PetscViewer for Exodus II input/output to use with the specified file
12966823f3c5SBlaise Bourdin 
12976823f3c5SBlaise Bourdin    Level: beginner
12986823f3c5SBlaise Bourdin 
12996823f3c5SBlaise Bourdin    Note:
13006823f3c5SBlaise Bourdin    This PetscViewer should be destroyed with PetscViewerDestroy().
13016823f3c5SBlaise Bourdin 
13026823f3c5SBlaise Bourdin .seealso: PetscViewerPushFormat(), PetscViewerDestroy(),
13036823f3c5SBlaise Bourdin           DMLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
13046823f3c5SBlaise Bourdin @*/
13056823f3c5SBlaise Bourdin PetscErrorCode PetscViewerExodusIIOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *exo)
13066823f3c5SBlaise Bourdin {
13076823f3c5SBlaise Bourdin   PetscErrorCode ierr;
13086823f3c5SBlaise Bourdin 
13096823f3c5SBlaise Bourdin   PetscFunctionBegin;
13106823f3c5SBlaise Bourdin   ierr = PetscViewerCreate(comm, exo);CHKERRQ(ierr);
13116823f3c5SBlaise Bourdin   ierr = PetscViewerSetType(*exo, PETSCVIEWEREXODUSII);CHKERRQ(ierr);
13126823f3c5SBlaise Bourdin   ierr = PetscViewerFileSetMode(*exo, type);CHKERRQ(ierr);
13136823f3c5SBlaise Bourdin   ierr = PetscViewerFileSetName(*exo, name);CHKERRQ(ierr);
13146823f3c5SBlaise Bourdin   ierr = PetscViewerSetFromOptions(*exo);CHKERRQ(ierr);
13151e50132fSMatthew G. Knepley   PetscFunctionReturn(0);
13161e50132fSMatthew G. Knepley }
13171e50132fSMatthew G. Knepley 
131833751fbdSMichael Lange /*@C
1319eaf898f9SPatrick Sanan   DMPlexCreateExodusFromFile - Create a DMPlex mesh from an ExodusII file.
132033751fbdSMichael Lange 
1321d083f849SBarry Smith   Collective
132233751fbdSMichael Lange 
132333751fbdSMichael Lange   Input Parameters:
132433751fbdSMichael Lange + comm  - The MPI communicator
132533751fbdSMichael Lange . filename - The name of the ExodusII file
132633751fbdSMichael Lange - interpolate - Create faces and edges in the mesh
132733751fbdSMichael Lange 
132833751fbdSMichael Lange   Output Parameter:
132933751fbdSMichael Lange . dm  - The DM object representing the mesh
133033751fbdSMichael Lange 
133133751fbdSMichael Lange   Level: beginner
133233751fbdSMichael Lange 
133333751fbdSMichael Lange .seealso: DMPLEX, DMCreate(), DMPlexCreateExodus()
133433751fbdSMichael Lange @*/
133533751fbdSMichael Lange PetscErrorCode DMPlexCreateExodusFromFile(MPI_Comm comm, const char filename[], PetscBool interpolate, DM *dm)
133633751fbdSMichael Lange {
133733751fbdSMichael Lange   PetscMPIInt    rank;
133833751fbdSMichael Lange   PetscErrorCode ierr;
133933751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
1340e84f7031SBlaise Bourdin   int   CPU_word_size = sizeof(PetscReal), IO_word_size = 0, exoid = -1;
134133751fbdSMichael Lange   float version;
134233751fbdSMichael Lange #endif
134333751fbdSMichael Lange 
134433751fbdSMichael Lange   PetscFunctionBegin;
134533751fbdSMichael Lange   PetscValidCharPointer(filename, 2);
1346ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
134733751fbdSMichael Lange #if defined(PETSC_HAVE_EXODUSII)
134833751fbdSMichael Lange   if (!rank) {
134933751fbdSMichael Lange     exoid = ex_open(filename, EX_READ, &CPU_word_size, &IO_word_size, &version);
135033751fbdSMichael Lange     if (exoid <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "ex_open(\"%s\",...) did not return a valid file ID", filename);
135133751fbdSMichael Lange   }
135233751fbdSMichael Lange   ierr = DMPlexCreateExodus(comm, exoid, interpolate, dm);CHKERRQ(ierr);
13536de834b4SMatthew G. Knepley   if (!rank) {PetscStackCallStandard(ex_close,(exoid));}
1354b458e8f1SJose E. Roman   PetscFunctionReturn(0);
135533751fbdSMichael Lange #else
135633751fbdSMichael Lange   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
135733751fbdSMichael Lange #endif
135833751fbdSMichael Lange }
135933751fbdSMichael Lange 
13608f861fbcSMatthew G. Knepley #if defined(PETSC_HAVE_EXODUSII)
13616823f3c5SBlaise Bourdin static PetscErrorCode ExodusGetCellType_Internal(const char *elem_type, DMPolytopeType *ct)
13628f861fbcSMatthew G. Knepley {
13638f861fbcSMatthew G. Knepley   PetscBool      flg;
13648f861fbcSMatthew G. Knepley   PetscErrorCode ierr;
13658f861fbcSMatthew G. Knepley 
13668f861fbcSMatthew G. Knepley   PetscFunctionBegin;
13678f861fbcSMatthew G. Knepley   *ct = DM_POLYTOPE_UNKNOWN;
13688f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "TRI", &flg);CHKERRQ(ierr);
13698f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
13708f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "TRI3", &flg);CHKERRQ(ierr);
13718f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRIANGLE; goto done;}
13728f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "QUAD", &flg);CHKERRQ(ierr);
13738f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13748f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "QUAD4", &flg);CHKERRQ(ierr);
13758f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13768f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "SHELL4", &flg);CHKERRQ(ierr);
13778f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_QUADRILATERAL; goto done;}
13788f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "TETRA", &flg);CHKERRQ(ierr);
13798f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
13808f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "TET4", &flg);CHKERRQ(ierr);
13818f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TETRAHEDRON; goto done;}
13828f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "WEDGE", &flg);CHKERRQ(ierr);
13838f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_TRI_PRISM; goto done;}
13848f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "HEX", &flg);CHKERRQ(ierr);
13858f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
13868f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "HEX8", &flg);CHKERRQ(ierr);
13878f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
13888f861fbcSMatthew G. Knepley   ierr = PetscStrcmp(elem_type, "HEXAHEDRON", &flg);CHKERRQ(ierr);
13898f861fbcSMatthew G. Knepley   if (flg) {*ct = DM_POLYTOPE_HEXAHEDRON; goto done;}
13908f861fbcSMatthew G. Knepley   SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Unrecognized element type %s", elem_type);
13918f861fbcSMatthew G. Knepley   done:
13928f861fbcSMatthew G. Knepley   PetscFunctionReturn(0);
13938f861fbcSMatthew G. Knepley }
13948f861fbcSMatthew G. Knepley #endif
13958f861fbcSMatthew G. Knepley 
1396552f7358SJed Brown /*@
139733751fbdSMichael Lange   DMPlexCreateExodus - Create a DMPlex mesh from an ExodusII file ID.
1398552f7358SJed Brown 
1399d083f849SBarry Smith   Collective
1400552f7358SJed Brown 
1401552f7358SJed Brown   Input Parameters:
1402552f7358SJed Brown + comm  - The MPI communicator
1403552f7358SJed Brown . exoid - The ExodusII id associated with a exodus file and obtained using ex_open
1404552f7358SJed Brown - interpolate - Create faces and edges in the mesh
1405552f7358SJed Brown 
1406552f7358SJed Brown   Output Parameter:
1407552f7358SJed Brown . dm  - The DM object representing the mesh
1408552f7358SJed Brown 
1409552f7358SJed Brown   Level: beginner
1410552f7358SJed Brown 
1411e84e7769SJed Brown .seealso: DMPLEX, DMCreate()
1412552f7358SJed Brown @*/
1413552f7358SJed Brown PetscErrorCode DMPlexCreateExodus(MPI_Comm comm, PetscInt exoid, PetscBool interpolate, DM *dm)
1414552f7358SJed Brown {
1415552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1416552f7358SJed Brown   PetscMPIInt    num_proc, rank;
1417ae9eebabSLisandro Dalcin   DMLabel        cellSets = NULL, faceSets = NULL, vertSets = NULL;
1418552f7358SJed Brown   PetscSection   coordSection;
1419552f7358SJed Brown   Vec            coordinates;
1420552f7358SJed Brown   PetscScalar    *coords;
1421552f7358SJed Brown   PetscInt       coordSize, v;
1422552f7358SJed Brown   PetscErrorCode ierr;
1423552f7358SJed Brown   /* Read from ex_get_init() */
1424552f7358SJed Brown   char title[PETSC_MAX_PATH_LEN+1];
14258f861fbcSMatthew G. Knepley   int  dim    = 0, dimEmbed = 0, numVertices = 0, numCells = 0, numHybridCells = 0;
1426552f7358SJed Brown   int  num_cs = 0, num_vs = 0, num_fs = 0;
1427552f7358SJed Brown #endif
1428552f7358SJed Brown 
1429552f7358SJed Brown   PetscFunctionBegin;
1430552f7358SJed Brown #if defined(PETSC_HAVE_EXODUSII)
1431ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm, &rank);CHKERRMPI(ierr);
1432ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm, &num_proc);CHKERRMPI(ierr);
1433552f7358SJed Brown   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
1434552f7358SJed Brown   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
1435552f7358SJed Brown   /* Open EXODUS II file and read basic informations on rank 0, then broadcast to all processors */
1436552f7358SJed Brown   if (!rank) {
1437580bdb30SBarry Smith     ierr = PetscMemzero(title,PETSC_MAX_PATH_LEN+1);CHKERRQ(ierr);
14388f861fbcSMatthew G. Knepley     PetscStackCallStandard(ex_get_init,(exoid, title, &dimEmbed, &numVertices, &numCells, &num_cs, &num_vs, &num_fs));
1439552f7358SJed Brown     if (!num_cs) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Exodus file does not contain any cell set\n");
1440552f7358SJed Brown   }
1441ffc4695bSBarry Smith   ierr = MPI_Bcast(title, PETSC_MAX_PATH_LEN+1, MPI_CHAR, 0, comm);CHKERRMPI(ierr);
1442ffc4695bSBarry Smith   ierr = MPI_Bcast(&dim, 1, MPI_INT, 0, comm);CHKERRMPI(ierr);
1443552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) *dm, title);CHKERRQ(ierr);
1444552f7358SJed Brown   ierr = DMPlexSetChart(*dm, 0, numCells+numVertices);CHKERRQ(ierr);
1445412e9a14SMatthew G. Knepley   /*   We do not want this label automatically computed, instead we compute it here */
1446412e9a14SMatthew G. Knepley   ierr = DMCreateLabel(*dm, "celltype");CHKERRQ(ierr);
1447552f7358SJed Brown 
1448552f7358SJed Brown   /* Read cell sets information */
1449552f7358SJed Brown   if (!rank) {
1450552f7358SJed Brown     PetscInt *cone;
1451e8f6893fSMatthew G. Knepley     int      c, cs, ncs, c_loc, v, v_loc;
1452552f7358SJed Brown     /* Read from ex_get_elem_blk_ids() */
1453e8f6893fSMatthew G. Knepley     int *cs_id, *cs_order;
1454552f7358SJed Brown     /* Read from ex_get_elem_block() */
1455552f7358SJed Brown     char buffer[PETSC_MAX_PATH_LEN+1];
1456e8f6893fSMatthew G. Knepley     int  num_cell_in_set, num_vertex_per_cell, num_hybrid, num_attr;
1457552f7358SJed Brown     /* Read from ex_get_elem_conn() */
1458552f7358SJed Brown     int *cs_connect;
1459552f7358SJed Brown 
1460552f7358SJed Brown     /* Get cell sets IDs */
1461e8f6893fSMatthew G. Knepley     ierr = PetscMalloc2(num_cs, &cs_id, num_cs, &cs_order);CHKERRQ(ierr);
14626de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_ELEM_BLOCK, cs_id));
1463552f7358SJed Brown     /* Read the cell set connectivity table and build mesh topology
1464552f7358SJed Brown        EXO standard requires that cells in cell sets be numbered sequentially and be pairwise disjoint. */
1465e8f6893fSMatthew G. Knepley     /* Check for a hybrid mesh */
1466e8f6893fSMatthew G. Knepley     for (cs = 0, num_hybrid = 0; cs < num_cs; ++cs) {
14678f861fbcSMatthew G. Knepley       DMPolytopeType ct;
14688f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
14698f861fbcSMatthew G. Knepley 
14708f861fbcSMatthew G. Knepley       ierr = PetscArrayzero(elem_type, sizeof(elem_type));CHKERRQ(ierr);
14718f861fbcSMatthew G. Knepley       PetscStackCallStandard(ex_get_elem_type,(exoid, cs_id[cs], elem_type));
14726823f3c5SBlaise Bourdin       ierr = ExodusGetCellType_Internal(elem_type, &ct);CHKERRQ(ierr);
14738f861fbcSMatthew G. Knepley       dim  = PetscMax(dim, DMPolytopeTypeGetDim(ct));
1474e8f6893fSMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
14758f861fbcSMatthew G. Knepley       switch (ct) {
14768f861fbcSMatthew G. Knepley         case DM_POLYTOPE_TRI_PRISM:
1477e8f6893fSMatthew G. Knepley           cs_order[cs] = cs;
1478e8f6893fSMatthew G. Knepley           numHybridCells += num_cell_in_set;
1479e8f6893fSMatthew G. Knepley           ++num_hybrid;
1480e8f6893fSMatthew G. Knepley           break;
1481e8f6893fSMatthew G. Knepley         default:
1482e8f6893fSMatthew G. Knepley           for (c = cs; c > cs-num_hybrid; --c) cs_order[c] = cs_order[c-1];
1483e8f6893fSMatthew G. Knepley           cs_order[cs-num_hybrid] = cs;
1484e8f6893fSMatthew G. Knepley       }
1485e8f6893fSMatthew G. Knepley     }
1486552f7358SJed Brown     /* First set sizes */
1487e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
14888f861fbcSMatthew G. Knepley       DMPolytopeType ct;
14898f861fbcSMatthew G. Knepley       char           elem_type[PETSC_MAX_PATH_LEN];
1490e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
14918f861fbcSMatthew G. Knepley 
14928f861fbcSMatthew G. Knepley       ierr = PetscArrayzero(elem_type, sizeof(elem_type));CHKERRQ(ierr);
14938f861fbcSMatthew G. Knepley       PetscStackCallStandard(ex_get_elem_type,(exoid, cs_id[cs], elem_type));
14946823f3c5SBlaise Bourdin       ierr = ExodusGetCellType_Internal(elem_type, &ct);CHKERRQ(ierr);
14956de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set,&num_vertex_per_cell, 0, 0, &num_attr));
1496552f7358SJed Brown       for (c_loc = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1497552f7358SJed Brown         ierr = DMPlexSetConeSize(*dm, c, num_vertex_per_cell);CHKERRQ(ierr);
1498412e9a14SMatthew G. Knepley         ierr = DMPlexSetCellType(*dm, c, ct);CHKERRQ(ierr);
1499552f7358SJed Brown       }
1500552f7358SJed Brown     }
1501412e9a14SMatthew G. Knepley     for (v = numCells; v < numCells+numVertices; ++v) {ierr = DMPlexSetCellType(*dm, v, DM_POLYTOPE_POINT);CHKERRQ(ierr);}
1502552f7358SJed Brown     ierr = DMSetUp(*dm);CHKERRQ(ierr);
1503e8f6893fSMatthew G. Knepley     for (ncs = 0, c = 0; ncs < num_cs; ++ncs) {
1504e8f6893fSMatthew G. Knepley       const PetscInt cs = cs_order[ncs];
15056de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_block,(exoid, EX_ELEM_BLOCK, cs_id[cs], buffer, &num_cell_in_set, &num_vertex_per_cell, 0, 0, &num_attr));
1506dcca6d9dSJed Brown       ierr = PetscMalloc2(num_vertex_per_cell*num_cell_in_set,&cs_connect,num_vertex_per_cell,&cone);CHKERRQ(ierr);
15076de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_conn,(exoid, EX_ELEM_BLOCK, cs_id[cs], cs_connect,NULL,NULL));
1508eaf898f9SPatrick Sanan       /* EXO uses Fortran-based indexing, DMPlex uses C-style and numbers cell first then vertices. */
1509552f7358SJed Brown       for (c_loc = 0, v = 0; c_loc < num_cell_in_set; ++c_loc, ++c) {
1510636e64ffSStefano Zampini         DMPolytopeType ct;
1511636e64ffSStefano Zampini 
1512552f7358SJed Brown         for (v_loc = 0; v_loc < num_vertex_per_cell; ++v_loc, ++v) {
1513552f7358SJed Brown           cone[v_loc] = cs_connect[v]+numCells-1;
1514552f7358SJed Brown         }
1515636e64ffSStefano Zampini         ierr = DMPlexGetCellType(*dm, c, &ct);CHKERRQ(ierr);
1516636e64ffSStefano Zampini         ierr = DMPlexInvertCell(ct, cone);CHKERRQ(ierr);
1517552f7358SJed Brown         ierr = DMPlexSetCone(*dm, c, cone);CHKERRQ(ierr);
1518ae9eebabSLisandro Dalcin         ierr = DMSetLabelValue_Fast(*dm, &cellSets, "Cell Sets", c, cs_id[cs]);CHKERRQ(ierr);
1519552f7358SJed Brown       }
1520552f7358SJed Brown       ierr = PetscFree2(cs_connect,cone);CHKERRQ(ierr);
1521552f7358SJed Brown     }
1522e8f6893fSMatthew G. Knepley     ierr = PetscFree2(cs_id, cs_order);CHKERRQ(ierr);
1523552f7358SJed Brown   }
15248f861fbcSMatthew G. Knepley   {
15258f861fbcSMatthew G. Knepley     PetscInt ints[] = {dim, dimEmbed};
15268f861fbcSMatthew G. Knepley 
1527ffc4695bSBarry Smith     ierr = MPI_Bcast(ints, 2, MPIU_INT, 0, comm);CHKERRMPI(ierr);
15288f861fbcSMatthew G. Knepley     ierr = DMSetDimension(*dm, ints[0]);CHKERRQ(ierr);
15298f861fbcSMatthew G. Knepley     ierr = DMSetCoordinateDim(*dm, ints[1]);CHKERRQ(ierr);
15308f861fbcSMatthew G. Knepley     dim      = ints[0];
15318f861fbcSMatthew G. Knepley     dimEmbed = ints[1];
15328f861fbcSMatthew G. Knepley   }
1533552f7358SJed Brown   ierr = DMPlexSymmetrize(*dm);CHKERRQ(ierr);
1534552f7358SJed Brown   ierr = DMPlexStratify(*dm);CHKERRQ(ierr);
1535552f7358SJed Brown   if (interpolate) {
15365fd9971aSMatthew G. Knepley     DM idm;
1537552f7358SJed Brown 
15389f074e33SMatthew G Knepley     ierr = DMPlexInterpolate(*dm, &idm);CHKERRQ(ierr);
1539552f7358SJed Brown     ierr = DMDestroy(dm);CHKERRQ(ierr);
1540552f7358SJed Brown     *dm  = idm;
1541552f7358SJed Brown   }
1542552f7358SJed Brown 
1543552f7358SJed Brown   /* Create vertex set label */
1544552f7358SJed Brown   if (!rank && (num_vs > 0)) {
1545552f7358SJed Brown     int vs, v;
1546552f7358SJed Brown     /* Read from ex_get_node_set_ids() */
1547552f7358SJed Brown     int *vs_id;
1548552f7358SJed Brown     /* Read from ex_get_node_set_param() */
1549f958083aSBlaise Bourdin     int num_vertex_in_set;
1550552f7358SJed Brown     /* Read from ex_get_node_set() */
1551552f7358SJed Brown     int *vs_vertex_list;
1552552f7358SJed Brown 
1553552f7358SJed Brown     /* Get vertex set ids */
1554785e854fSJed Brown     ierr = PetscMalloc1(num_vs, &vs_id);CHKERRQ(ierr);
15556de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_NODE_SET, vs_id));
1556552f7358SJed Brown     for (vs = 0; vs < num_vs; ++vs) {
15576de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_NODE_SET, vs_id[vs], &num_vertex_in_set, NULL));
1558785e854fSJed Brown       ierr = PetscMalloc1(num_vertex_in_set, &vs_vertex_list);CHKERRQ(ierr);
15596de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set,(exoid, EX_NODE_SET, vs_id[vs], vs_vertex_list, NULL));
1560552f7358SJed Brown       for (v = 0; v < num_vertex_in_set; ++v) {
1561ae9eebabSLisandro Dalcin         ierr = DMSetLabelValue_Fast(*dm, &vertSets, "Vertex Sets", vs_vertex_list[v]+numCells-1, vs_id[vs]);CHKERRQ(ierr);
1562552f7358SJed Brown       }
1563552f7358SJed Brown       ierr = PetscFree(vs_vertex_list);CHKERRQ(ierr);
1564552f7358SJed Brown     }
1565552f7358SJed Brown     ierr = PetscFree(vs_id);CHKERRQ(ierr);
1566552f7358SJed Brown   }
1567552f7358SJed Brown   /* Read coordinates */
156869d8a9ceSMatthew G. Knepley   ierr = DMGetCoordinateSection(*dm, &coordSection);CHKERRQ(ierr);
1569552f7358SJed Brown   ierr = PetscSectionSetNumFields(coordSection, 1);CHKERRQ(ierr);
15708f861fbcSMatthew G. Knepley   ierr = PetscSectionSetFieldComponents(coordSection, 0, dimEmbed);CHKERRQ(ierr);
1571552f7358SJed Brown   ierr = PetscSectionSetChart(coordSection, numCells, numCells + numVertices);CHKERRQ(ierr);
1572552f7358SJed Brown   for (v = numCells; v < numCells+numVertices; ++v) {
15738f861fbcSMatthew G. Knepley     ierr = PetscSectionSetDof(coordSection, v, dimEmbed);CHKERRQ(ierr);
15748f861fbcSMatthew G. Knepley     ierr = PetscSectionSetFieldDof(coordSection, v, 0, dimEmbed);CHKERRQ(ierr);
1575552f7358SJed Brown   }
1576552f7358SJed Brown   ierr = PetscSectionSetUp(coordSection);CHKERRQ(ierr);
1577552f7358SJed Brown   ierr = PetscSectionGetStorageSize(coordSection, &coordSize);CHKERRQ(ierr);
15788b9ced59SLisandro Dalcin   ierr = VecCreate(PETSC_COMM_SELF, &coordinates);CHKERRQ(ierr);
1579552f7358SJed Brown   ierr = PetscObjectSetName((PetscObject) coordinates, "coordinates");CHKERRQ(ierr);
1580552f7358SJed Brown   ierr = VecSetSizes(coordinates, coordSize, PETSC_DETERMINE);CHKERRQ(ierr);
15818f861fbcSMatthew G. Knepley   ierr = VecSetBlockSize(coordinates, dimEmbed);CHKERRQ(ierr);
15822eb5907fSJed Brown   ierr = VecSetType(coordinates,VECSTANDARD);CHKERRQ(ierr);
1583552f7358SJed Brown   ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr);
15840aba08f6SKarl Rupp   if (!rank) {
1585e84f7031SBlaise Bourdin     PetscReal *x, *y, *z;
1586552f7358SJed Brown 
1587dcca6d9dSJed Brown     ierr = PetscMalloc3(numVertices,&x,numVertices,&y,numVertices,&z);CHKERRQ(ierr);
15886de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_coord,(exoid, x, y, z));
15898f861fbcSMatthew G. Knepley     if (dimEmbed > 0) {
15908f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+0] = x[v];
15910d644c17SKarl Rupp     }
15928f861fbcSMatthew G. Knepley     if (dimEmbed > 1) {
15938f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+1] = y[v];
15940d644c17SKarl Rupp     }
15958f861fbcSMatthew G. Knepley     if (dimEmbed > 2) {
15968f861fbcSMatthew G. Knepley       for (v = 0; v < numVertices; ++v) coords[v*dimEmbed+2] = z[v];
15970d644c17SKarl Rupp     }
1598552f7358SJed Brown     ierr = PetscFree3(x,y,z);CHKERRQ(ierr);
1599552f7358SJed Brown   }
1600552f7358SJed Brown   ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr);
1601552f7358SJed Brown   ierr = DMSetCoordinatesLocal(*dm, coordinates);CHKERRQ(ierr);
1602552f7358SJed Brown   ierr = VecDestroy(&coordinates);CHKERRQ(ierr);
1603552f7358SJed Brown 
1604552f7358SJed Brown   /* Create side set label */
1605552f7358SJed Brown   if (!rank && interpolate && (num_fs > 0)) {
1606552f7358SJed Brown     int fs, f, voff;
1607552f7358SJed Brown     /* Read from ex_get_side_set_ids() */
1608552f7358SJed Brown     int *fs_id;
1609552f7358SJed Brown     /* Read from ex_get_side_set_param() */
1610f958083aSBlaise Bourdin     int num_side_in_set;
1611552f7358SJed Brown     /* Read from ex_get_side_set_node_list() */
1612552f7358SJed Brown     int *fs_vertex_count_list, *fs_vertex_list;
1613ef073283Smichael_afanasiev     /* Read side set labels */
1614ef073283Smichael_afanasiev     char fs_name[MAX_STR_LENGTH+1];
1615552f7358SJed Brown 
1616552f7358SJed Brown     /* Get side set ids */
1617785e854fSJed Brown     ierr = PetscMalloc1(num_fs, &fs_id);CHKERRQ(ierr);
16186de834b4SMatthew G. Knepley     PetscStackCallStandard(ex_get_ids,(exoid, EX_SIDE_SET, fs_id));
1619552f7358SJed Brown     for (fs = 0; fs < num_fs; ++fs) {
16206de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_set_param,(exoid, EX_SIDE_SET, fs_id[fs], &num_side_in_set, NULL));
1621dcca6d9dSJed Brown       ierr = PetscMalloc2(num_side_in_set,&fs_vertex_count_list,num_side_in_set*4,&fs_vertex_list);CHKERRQ(ierr);
16226de834b4SMatthew G. Knepley       PetscStackCallStandard(ex_get_side_set_node_list,(exoid, fs_id[fs], fs_vertex_count_list, fs_vertex_list));
1623ef073283Smichael_afanasiev       /* Get the specific name associated with this side set ID. */
1624ef073283Smichael_afanasiev       int fs_name_err = ex_get_name(exoid, EX_SIDE_SET, fs_id[fs], fs_name);
1625552f7358SJed Brown       for (f = 0, voff = 0; f < num_side_in_set; ++f) {
16260298fd71SBarry Smith         const PetscInt *faces   = NULL;
1627552f7358SJed Brown         PetscInt       faceSize = fs_vertex_count_list[f], numFaces;
1628552f7358SJed Brown         PetscInt       faceVertices[4], v;
1629552f7358SJed Brown 
1630552f7358SJed Brown         if (faceSize > 4) SETERRQ1(comm, PETSC_ERR_ARG_WRONG, "ExodusII side cannot have %d > 4 vertices", faceSize);
1631552f7358SJed Brown         for (v = 0; v < faceSize; ++v, ++voff) {
1632552f7358SJed Brown           faceVertices[v] = fs_vertex_list[voff]+numCells-1;
1633552f7358SJed Brown         }
1634552f7358SJed Brown         ierr = DMPlexGetFullJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1635552f7358SJed Brown         if (numFaces != 1) SETERRQ3(comm, PETSC_ERR_ARG_WRONG, "Invalid ExodusII side %d in set %d maps to %d faces", f, fs, numFaces);
1636ae9eebabSLisandro Dalcin         ierr = DMSetLabelValue_Fast(*dm, &faceSets, "Face Sets", faces[0], fs_id[fs]);CHKERRQ(ierr);
1637ef073283Smichael_afanasiev         /* Only add the label if one has been detected for this side set. */
1638ef073283Smichael_afanasiev         if (!fs_name_err) {
1639ef073283Smichael_afanasiev           ierr = DMSetLabelValue(*dm, fs_name, faces[0], fs_id[fs]);CHKERRQ(ierr);
1640ef073283Smichael_afanasiev         }
1641552f7358SJed Brown         ierr = DMPlexRestoreJoin(*dm, faceSize, faceVertices, &numFaces, &faces);CHKERRQ(ierr);
1642552f7358SJed Brown       }
1643552f7358SJed Brown       ierr = PetscFree2(fs_vertex_count_list,fs_vertex_list);CHKERRQ(ierr);
1644552f7358SJed Brown     }
1645552f7358SJed Brown     ierr = PetscFree(fs_id);CHKERRQ(ierr);
1646552f7358SJed Brown   }
1647ae9eebabSLisandro Dalcin 
1648ae9eebabSLisandro Dalcin   { /* Create Cell/Face/Vertex Sets labels at all processes */
1649ae9eebabSLisandro Dalcin     enum {n = 3};
1650ae9eebabSLisandro Dalcin     PetscBool flag[n];
1651ae9eebabSLisandro Dalcin 
1652ae9eebabSLisandro Dalcin     flag[0] = cellSets ? PETSC_TRUE : PETSC_FALSE;
1653ae9eebabSLisandro Dalcin     flag[1] = faceSets ? PETSC_TRUE : PETSC_FALSE;
1654ae9eebabSLisandro Dalcin     flag[2] = vertSets ? PETSC_TRUE : PETSC_FALSE;
1655ae9eebabSLisandro Dalcin     ierr = MPI_Bcast(flag, n, MPIU_BOOL, 0, comm);CHKERRMPI(ierr);
1656ae9eebabSLisandro Dalcin     if (flag[0]) {ierr = DMCreateLabel(*dm, "Cell Sets");CHKERRQ(ierr);}
1657ae9eebabSLisandro Dalcin     if (flag[1]) {ierr = DMCreateLabel(*dm, "Face Sets");CHKERRQ(ierr);}
1658ae9eebabSLisandro Dalcin     if (flag[2]) {ierr = DMCreateLabel(*dm, "Vertex Sets");CHKERRQ(ierr);}
1659ae9eebabSLisandro Dalcin   }
1660b458e8f1SJose E. Roman   PetscFunctionReturn(0);
1661552f7358SJed Brown #else
1662552f7358SJed Brown   SETERRQ(comm, PETSC_ERR_SUP, "This method requires ExodusII support. Reconfigure using --download-exodusii");
1663552f7358SJed Brown #endif
1664552f7358SJed Brown }
1665