xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 729ab6d9824aff21597cdd01187d6d9cd1cdf05e)
15c6c1daeSBarry Smith #include <petsc-private/viewerimpl.h>    /*I   "petscsys.h"   I*/
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
35c6c1daeSBarry Smith 
45c6c1daeSBarry Smith typedef struct GroupList {
55c6c1daeSBarry Smith   const char       *name;
65c6c1daeSBarry Smith   struct GroupList *next;
75c6c1daeSBarry Smith } GroupList;
85c6c1daeSBarry Smith 
95c6c1daeSBarry Smith typedef struct {
105c6c1daeSBarry Smith   char          *filename;
115c6c1daeSBarry Smith   PetscFileMode btype;
125c6c1daeSBarry Smith   hid_t         file_id;
135c6c1daeSBarry Smith   PetscInt      timestep;
145c6c1daeSBarry Smith   GroupList     *groups;
155c6c1daeSBarry Smith } PetscViewer_HDF5;
165c6c1daeSBarry Smith 
175c6c1daeSBarry Smith #undef __FUNCT__
185c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileClose_HDF5"
195c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
205c6c1daeSBarry Smith {
215c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
225c6c1daeSBarry Smith   PetscErrorCode   ierr;
235c6c1daeSBarry Smith 
245c6c1daeSBarry Smith   PetscFunctionBegin;
255c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
26*729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
275c6c1daeSBarry Smith   PetscFunctionReturn(0);
285c6c1daeSBarry Smith }
295c6c1daeSBarry Smith 
305c6c1daeSBarry Smith #undef __FUNCT__
315c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerDestroy_HDF5"
325c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
335c6c1daeSBarry Smith {
345c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
355c6c1daeSBarry Smith   PetscErrorCode   ierr;
365c6c1daeSBarry Smith 
375c6c1daeSBarry Smith   PetscFunctionBegin;
385c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
395c6c1daeSBarry Smith   while (hdf5->groups) {
405c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
415c6c1daeSBarry Smith 
425c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
435c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
445c6c1daeSBarry Smith     hdf5->groups = tmp;
455c6c1daeSBarry Smith   }
465c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
470b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
480b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
495c6c1daeSBarry Smith   PetscFunctionReturn(0);
505c6c1daeSBarry Smith }
515c6c1daeSBarry Smith 
525c6c1daeSBarry Smith #undef __FUNCT__
535c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetMode_HDF5"
545c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
555c6c1daeSBarry Smith {
565c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
575c6c1daeSBarry Smith 
585c6c1daeSBarry Smith   PetscFunctionBegin;
595c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
605c6c1daeSBarry Smith   hdf5->btype = type;
615c6c1daeSBarry Smith   PetscFunctionReturn(0);
625c6c1daeSBarry Smith }
635c6c1daeSBarry Smith 
645c6c1daeSBarry Smith #undef __FUNCT__
655c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerFileSetName_HDF5"
665c6c1daeSBarry Smith PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
675c6c1daeSBarry Smith {
685c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
695c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
705c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
715c6c1daeSBarry Smith #endif
725c6c1daeSBarry Smith   hid_t             plist_id;
735c6c1daeSBarry Smith   PetscErrorCode    ierr;
745c6c1daeSBarry Smith 
755c6c1daeSBarry Smith   PetscFunctionBegin;
765c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
775c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
78*729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
795c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
80*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
815c6c1daeSBarry Smith #endif
825c6c1daeSBarry Smith   /* Create or open the file collectively */
835c6c1daeSBarry Smith   switch (hdf5->btype) {
845c6c1daeSBarry Smith   case FILE_MODE_READ:
85*729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
865c6c1daeSBarry Smith     break;
875c6c1daeSBarry Smith   case FILE_MODE_APPEND:
88*729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
895c6c1daeSBarry Smith     break;
905c6c1daeSBarry Smith   case FILE_MODE_WRITE:
91*729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
925c6c1daeSBarry Smith     break;
935c6c1daeSBarry Smith   default:
945c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
955c6c1daeSBarry Smith   }
965c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
97*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
985c6c1daeSBarry Smith   PetscFunctionReturn(0);
995c6c1daeSBarry Smith }
1005c6c1daeSBarry Smith 
1015c6c1daeSBarry Smith #undef __FUNCT__
1025c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerCreate_HDF5"
1038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
1045c6c1daeSBarry Smith {
1055c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
1065c6c1daeSBarry Smith   PetscErrorCode   ierr;
1075c6c1daeSBarry Smith 
1085c6c1daeSBarry Smith   PetscFunctionBegin;
109b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
1105c6c1daeSBarry Smith 
1115c6c1daeSBarry Smith   v->data         = (void*) hdf5;
1125c6c1daeSBarry Smith   v->ops->destroy = PetscViewerDestroy_HDF5;
1135c6c1daeSBarry Smith   v->ops->flush   = 0;
1145c6c1daeSBarry Smith   hdf5->btype     = (PetscFileMode) -1;
1155c6c1daeSBarry Smith   hdf5->filename  = 0;
1165c6c1daeSBarry Smith   hdf5->timestep  = -1;
1170298fd71SBarry Smith   hdf5->groups    = NULL;
1185c6c1daeSBarry Smith 
1190b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
1200b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
1215c6c1daeSBarry Smith   PetscFunctionReturn(0);
1225c6c1daeSBarry Smith }
1235c6c1daeSBarry Smith 
1245c6c1daeSBarry Smith #undef __FUNCT__
1255c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5Open"
1265c6c1daeSBarry Smith /*@C
1275c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
1285c6c1daeSBarry Smith 
1295c6c1daeSBarry Smith    Collective on MPI_Comm
1305c6c1daeSBarry Smith 
1315c6c1daeSBarry Smith    Input Parameters:
1325c6c1daeSBarry Smith +  comm - MPI communicator
1335c6c1daeSBarry Smith .  name - name of file
1345c6c1daeSBarry Smith -  type - type of file
1355c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
1365c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
1375c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
1385c6c1daeSBarry Smith 
1395c6c1daeSBarry Smith    Output Parameter:
1405c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
1415c6c1daeSBarry Smith 
1425c6c1daeSBarry Smith    Level: beginner
1435c6c1daeSBarry Smith 
1445c6c1daeSBarry Smith    Note:
1455c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
1465c6c1daeSBarry Smith 
1475c6c1daeSBarry Smith    Concepts: HDF5 files
1485c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
1495c6c1daeSBarry Smith 
1505c6c1daeSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
1515c6c1daeSBarry Smith           VecView(), MatView(), VecLoad(), MatLoad(),
1525c6c1daeSBarry Smith           PetscFileMode, PetscViewer
1535c6c1daeSBarry Smith @*/
1545c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
1555c6c1daeSBarry Smith {
1565c6c1daeSBarry Smith   PetscErrorCode ierr;
1575c6c1daeSBarry Smith 
1585c6c1daeSBarry Smith   PetscFunctionBegin;
1595c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
1605c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
1615c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
1625c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
1635c6c1daeSBarry Smith   PetscFunctionReturn(0);
1645c6c1daeSBarry Smith }
1655c6c1daeSBarry Smith 
1665c6c1daeSBarry Smith #undef __FUNCT__
1675c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetFileId"
1685c6c1daeSBarry Smith /*@C
1695c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
1705c6c1daeSBarry Smith 
1715c6c1daeSBarry Smith   Not collective
1725c6c1daeSBarry Smith 
1735c6c1daeSBarry Smith   Input Parameter:
1745c6c1daeSBarry Smith . viewer - the PetscViewer
1755c6c1daeSBarry Smith 
1765c6c1daeSBarry Smith   Output Parameter:
1775c6c1daeSBarry Smith . file_id - The file id
1785c6c1daeSBarry Smith 
1795c6c1daeSBarry Smith   Level: intermediate
1805c6c1daeSBarry Smith 
1815c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
1825c6c1daeSBarry Smith @*/
1835c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
1845c6c1daeSBarry Smith {
1855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1865c6c1daeSBarry Smith 
1875c6c1daeSBarry Smith   PetscFunctionBegin;
1885c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1895c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
1905c6c1daeSBarry Smith   PetscFunctionReturn(0);
1915c6c1daeSBarry Smith }
1925c6c1daeSBarry Smith 
1935c6c1daeSBarry Smith #undef __FUNCT__
1945c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PushGroup"
1955c6c1daeSBarry Smith /*@C
1965c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
1975c6c1daeSBarry Smith 
1985c6c1daeSBarry Smith   Not collective
1995c6c1daeSBarry Smith 
2005c6c1daeSBarry Smith   Input Parameters:
2015c6c1daeSBarry Smith + viewer - the PetscViewer
2025c6c1daeSBarry Smith - name - The group name
2035c6c1daeSBarry Smith 
2045c6c1daeSBarry Smith   Level: intermediate
2055c6c1daeSBarry Smith 
2065c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
2075c6c1daeSBarry Smith @*/
2085c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
2095c6c1daeSBarry Smith {
2105c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2115c6c1daeSBarry Smith   GroupList        *groupNode;
2125c6c1daeSBarry Smith   PetscErrorCode   ierr;
2135c6c1daeSBarry Smith 
2145c6c1daeSBarry Smith   PetscFunctionBegin;
2155c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2165c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
2175c6c1daeSBarry Smith   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
2185c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
219a297a907SKarl Rupp 
2205c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
2215c6c1daeSBarry Smith   hdf5->groups    = groupNode;
2225c6c1daeSBarry Smith   PetscFunctionReturn(0);
2235c6c1daeSBarry Smith }
2245c6c1daeSBarry Smith 
2255c6c1daeSBarry Smith #undef __FUNCT__
2265c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5PopGroup"
2273ef9c667SSatish Balay /*@
2285c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
2295c6c1daeSBarry Smith 
2305c6c1daeSBarry Smith   Not collective
2315c6c1daeSBarry Smith 
2325c6c1daeSBarry Smith   Input Parameter:
2335c6c1daeSBarry Smith . viewer - the PetscViewer
2345c6c1daeSBarry Smith 
2355c6c1daeSBarry Smith   Level: intermediate
2365c6c1daeSBarry Smith 
2375c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
2385c6c1daeSBarry Smith @*/
2395c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
2405c6c1daeSBarry Smith {
2415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2425c6c1daeSBarry Smith   GroupList        *groupNode;
2435c6c1daeSBarry Smith   PetscErrorCode   ierr;
2445c6c1daeSBarry Smith 
2455c6c1daeSBarry Smith   PetscFunctionBegin;
2465c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
24782f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
2485c6c1daeSBarry Smith   groupNode    = hdf5->groups;
2495c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
2505c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
2515c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
2525c6c1daeSBarry Smith   PetscFunctionReturn(0);
2535c6c1daeSBarry Smith }
2545c6c1daeSBarry Smith 
2555c6c1daeSBarry Smith #undef __FUNCT__
2565c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetGroup"
2575c6c1daeSBarry Smith /*@C
2580298fd71SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
2595c6c1daeSBarry Smith 
2605c6c1daeSBarry Smith   Not collective
2615c6c1daeSBarry Smith 
2625c6c1daeSBarry Smith   Input Parameter:
2635c6c1daeSBarry Smith . viewer - the PetscViewer
2645c6c1daeSBarry Smith 
2655c6c1daeSBarry Smith   Output Parameter:
2665c6c1daeSBarry Smith . name - The group name
2675c6c1daeSBarry Smith 
2685c6c1daeSBarry Smith   Level: intermediate
2695c6c1daeSBarry Smith 
2705c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
2715c6c1daeSBarry Smith @*/
2725c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
2735c6c1daeSBarry Smith {
2745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
2755c6c1daeSBarry Smith 
2765c6c1daeSBarry Smith   PetscFunctionBegin;
2775c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
278c959eef4SJed Brown   PetscValidPointer(name,2);
279a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
2800298fd71SBarry Smith   else *name = NULL;
2815c6c1daeSBarry Smith   PetscFunctionReturn(0);
2825c6c1daeSBarry Smith }
2835c6c1daeSBarry Smith 
2845c6c1daeSBarry Smith #undef __FUNCT__
2855c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
2863ef9c667SSatish Balay /*@
2875c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
2885c6c1daeSBarry Smith 
2895c6c1daeSBarry Smith   Not collective
2905c6c1daeSBarry Smith 
2915c6c1daeSBarry Smith   Input Parameter:
2925c6c1daeSBarry Smith . viewer - the PetscViewer
2935c6c1daeSBarry Smith 
2945c6c1daeSBarry Smith   Level: intermediate
2955c6c1daeSBarry Smith 
2965c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
2975c6c1daeSBarry Smith @*/
2985c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
2995c6c1daeSBarry Smith {
3005c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3015c6c1daeSBarry Smith 
3025c6c1daeSBarry Smith   PetscFunctionBegin;
3035c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3045c6c1daeSBarry Smith   ++hdf5->timestep;
3055c6c1daeSBarry Smith   PetscFunctionReturn(0);
3065c6c1daeSBarry Smith }
3075c6c1daeSBarry Smith 
3085c6c1daeSBarry Smith #undef __FUNCT__
3095c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5SetTimestep"
3103ef9c667SSatish Balay /*@
3115c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
3125c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
3135c6c1daeSBarry Smith 
3145c6c1daeSBarry Smith   Not collective
3155c6c1daeSBarry Smith 
3165c6c1daeSBarry Smith   Input Parameters:
3175c6c1daeSBarry Smith + viewer - the PetscViewer
3185c6c1daeSBarry Smith - timestep - The timestep number
3195c6c1daeSBarry Smith 
3205c6c1daeSBarry Smith   Level: intermediate
3215c6c1daeSBarry Smith 
3225c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
3235c6c1daeSBarry Smith @*/
3245c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
3255c6c1daeSBarry Smith {
3265c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3275c6c1daeSBarry Smith 
3285c6c1daeSBarry Smith   PetscFunctionBegin;
3295c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3305c6c1daeSBarry Smith   hdf5->timestep = timestep;
3315c6c1daeSBarry Smith   PetscFunctionReturn(0);
3325c6c1daeSBarry Smith }
3335c6c1daeSBarry Smith 
3345c6c1daeSBarry Smith #undef __FUNCT__
3355c6c1daeSBarry Smith #define __FUNCT__ "PetscViewerHDF5GetTimestep"
3363ef9c667SSatish Balay /*@
3375c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
3385c6c1daeSBarry Smith 
3395c6c1daeSBarry Smith   Not collective
3405c6c1daeSBarry Smith 
3415c6c1daeSBarry Smith   Input Parameter:
3425c6c1daeSBarry Smith . viewer - the PetscViewer
3435c6c1daeSBarry Smith 
3445c6c1daeSBarry Smith   Output Parameter:
3455c6c1daeSBarry Smith . timestep - The timestep number
3465c6c1daeSBarry Smith 
3475c6c1daeSBarry Smith   Level: intermediate
3485c6c1daeSBarry Smith 
3495c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
3505c6c1daeSBarry Smith @*/
3515c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
3525c6c1daeSBarry Smith {
3535c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3545c6c1daeSBarry Smith 
3555c6c1daeSBarry Smith   PetscFunctionBegin;
3565c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
3575c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
3585c6c1daeSBarry Smith   *timestep = hdf5->timestep;
3595c6c1daeSBarry Smith   PetscFunctionReturn(0);
3605c6c1daeSBarry Smith }
3615c6c1daeSBarry Smith 
36236bce6e8SMatthew G. Knepley #undef __FUNCT__
36336bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscDataTypeToHDF5DataType"
36436bce6e8SMatthew G. Knepley /*@C
36536bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
36636bce6e8SMatthew G. Knepley 
36736bce6e8SMatthew G. Knepley   Not collective
36836bce6e8SMatthew G. Knepley 
36936bce6e8SMatthew G. Knepley   Input Parameter:
37036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
37136bce6e8SMatthew G. Knepley 
37236bce6e8SMatthew G. Knepley   Output Parameter:
37336bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
37436bce6e8SMatthew G. Knepley 
37536bce6e8SMatthew G. Knepley   Level: advanced
37636bce6e8SMatthew G. Knepley 
37736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
37836bce6e8SMatthew G. Knepley @*/
37936bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
38036bce6e8SMatthew G. Knepley {
38136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
38236bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
38336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
38436bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
38536bce6e8SMatthew G. Knepley #else
38636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
38736bce6e8SMatthew G. Knepley #endif
38836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
38936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
39036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
39136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
39236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
39336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
39436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
39536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
3967e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
39736bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
39836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
39936bce6e8SMatthew G. Knepley }
40036bce6e8SMatthew G. Knepley 
40136bce6e8SMatthew G. Knepley #undef __FUNCT__
40236bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
40336bce6e8SMatthew G. Knepley /*@C
40436bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
40536bce6e8SMatthew G. Knepley 
40636bce6e8SMatthew G. Knepley   Not collective
40736bce6e8SMatthew G. Knepley 
40836bce6e8SMatthew G. Knepley   Input Parameter:
40936bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
41036bce6e8SMatthew G. Knepley 
41136bce6e8SMatthew G. Knepley   Output Parameter:
41236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
41336bce6e8SMatthew G. Knepley 
41436bce6e8SMatthew G. Knepley   Level: advanced
41536bce6e8SMatthew G. Knepley 
41636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
41736bce6e8SMatthew G. Knepley @*/
41836bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
41936bce6e8SMatthew G. Knepley {
42036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
42136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
42236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
42336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
42436bce6e8SMatthew G. Knepley #else
42536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
42636bce6e8SMatthew G. Knepley #endif
42736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
42836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
42936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
43036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
43136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
43236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
4337e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
43436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
43536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
43636bce6e8SMatthew G. Knepley }
43736bce6e8SMatthew G. Knepley 
43836bce6e8SMatthew G. Knepley #undef __FUNCT__
43936bce6e8SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
44036bce6e8SMatthew G. Knepley /*@
44136bce6e8SMatthew G. Knepley  PetscViewerHDF5WriteAttribute - Write a scalar attribute
44236bce6e8SMatthew G. Knepley 
44336bce6e8SMatthew G. Knepley   Input Parameters:
44436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
44536bce6e8SMatthew G. Knepley . parent - The parent name
44636bce6e8SMatthew G. Knepley . name   - The attribute name
44736bce6e8SMatthew G. Knepley . datatype - The attribute type
44836bce6e8SMatthew G. Knepley - value    - The attribute value
44936bce6e8SMatthew G. Knepley 
45036bce6e8SMatthew G. Knepley   Level: advanced
45136bce6e8SMatthew G. Knepley 
452e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
45336bce6e8SMatthew G. Knepley @*/
45436bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
45536bce6e8SMatthew G. Knepley {
456*729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
45736bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
45836bce6e8SMatthew G. Knepley 
45936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
4605cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
46136bce6e8SMatthew G. Knepley   PetscValidPointer(parent, 2);
46236bce6e8SMatthew G. Knepley   PetscValidPointer(name, 3);
46336bce6e8SMatthew G. Knepley   PetscValidPointer(value, 4);
46436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
4657e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
4667e97c476SMatthew G. Knepley     size_t len;
4677e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
468*729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
4697e97c476SMatthew G. Knepley   }
47036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
471*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
47236bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
473*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
474*729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
47536bce6e8SMatthew G. Knepley #else
476*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
477*729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
47836bce6e8SMatthew G. Knepley #endif
479*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
480*729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
481*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
482*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
483*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
48436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
48536bce6e8SMatthew G. Knepley }
48636bce6e8SMatthew G. Knepley 
487ad0c61feSMatthew G. Knepley #undef __FUNCT__
488ad0c61feSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5ReadAttribute"
489ad0c61feSMatthew G. Knepley /*@
490ad0c61feSMatthew G. Knepley  PetscViewerHDF5ReadAttribute - Read a scalar attribute
491ad0c61feSMatthew G. Knepley 
492ad0c61feSMatthew G. Knepley   Input Parameters:
493ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
494ad0c61feSMatthew G. Knepley . parent - The parent name
495ad0c61feSMatthew G. Knepley . name   - The attribute name
496ad0c61feSMatthew G. Knepley - datatype - The attribute type
497ad0c61feSMatthew G. Knepley 
498ad0c61feSMatthew G. Knepley   Output Parameter:
499ad0c61feSMatthew G. Knepley . value    - The attribute value
500ad0c61feSMatthew G. Knepley 
501ad0c61feSMatthew G. Knepley   Level: advanced
502ad0c61feSMatthew G. Knepley 
503e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
504ad0c61feSMatthew G. Knepley @*/
505ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
506ad0c61feSMatthew G. Knepley {
507*729ab6d9SBarry Smith   hid_t          h5, dataspace, dataset, attribute, dtype;
508ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
509ad0c61feSMatthew G. Knepley 
510ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
5115cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
512ad0c61feSMatthew G. Knepley   PetscValidPointer(parent, 2);
513ad0c61feSMatthew G. Knepley   PetscValidPointer(name, 3);
514ad0c61feSMatthew G. Knepley   PetscValidPointer(value, 4);
515ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
516ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
517*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
518ad0c61feSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
519*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
520ad0c61feSMatthew G. Knepley #else
521*729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
522ad0c61feSMatthew G. Knepley #endif
523*729ab6d9SBarry Smith   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
52470efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
525*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
526*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Dclose,(dataset));
527*729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
528ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
529ad0c61feSMatthew G. Knepley }
530ad0c61feSMatthew G. Knepley 
531e7bdbf8eSMatthew G. Knepley #undef __FUNCT__
5325cdeb108SMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasObject"
5335cdeb108SMatthew G. Knepley static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
5345cdeb108SMatthew G. Knepley {
5355cdeb108SMatthew G. Knepley   hid_t          h5;
5365cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
5375cdeb108SMatthew G. Knepley 
5385cdeb108SMatthew G. Knepley   PetscFunctionBegin;
5395cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5405cdeb108SMatthew G. Knepley   PetscValidPointer(name, 2);
5415cdeb108SMatthew G. Knepley   PetscValidPointer(has, 3);
5425cdeb108SMatthew G. Knepley   *has = PETSC_FALSE;
543c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
5445cdeb108SMatthew G. Knepley   if (H5Lexists(h5, name, H5P_DEFAULT)) {
5455cdeb108SMatthew G. Knepley     H5O_info_t info;
5465cdeb108SMatthew G. Knepley     hid_t      obj;
5475cdeb108SMatthew G. Knepley 
548*729ab6d9SBarry Smith     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
549*729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oget_info,(obj, &info));
5505cdeb108SMatthew G. Knepley     if (otype == info.type) *has = PETSC_TRUE;
551*729ab6d9SBarry Smith     PetscStackCallHDF5(H5Oclose,(obj));
5525cdeb108SMatthew G. Knepley   }
5535cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
5545cdeb108SMatthew G. Knepley }
5555cdeb108SMatthew G. Knepley 
5565cdeb108SMatthew G. Knepley #undef __FUNCT__
557e7bdbf8eSMatthew G. Knepley #define __FUNCT__ "PetscViewerHDF5HasAttribute"
558e7bdbf8eSMatthew G. Knepley /*@
559e7bdbf8eSMatthew G. Knepley  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
560e7bdbf8eSMatthew G. Knepley 
561e7bdbf8eSMatthew G. Knepley   Input Parameters:
562e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
563e7bdbf8eSMatthew G. Knepley . parent - The parent name
564e7bdbf8eSMatthew G. Knepley - name   - The attribute name
565e7bdbf8eSMatthew G. Knepley 
566e7bdbf8eSMatthew G. Knepley   Output Parameter:
567e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
568e7bdbf8eSMatthew G. Knepley 
569e7bdbf8eSMatthew G. Knepley   Level: advanced
570e7bdbf8eSMatthew G. Knepley 
571e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
572e7bdbf8eSMatthew G. Knepley @*/
573e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
574e7bdbf8eSMatthew G. Knepley {
575*729ab6d9SBarry Smith   hid_t          h5, dataset;
5765cdeb108SMatthew G. Knepley   htri_t         hhas;
5775cdeb108SMatthew G. Knepley   PetscBool      exists;
578e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
579e7bdbf8eSMatthew G. Knepley 
580e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
5815cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
582e7bdbf8eSMatthew G. Knepley   PetscValidPointer(parent, 2);
583e7bdbf8eSMatthew G. Knepley   PetscValidPointer(name, 3);
584e7bdbf8eSMatthew G. Knepley   PetscValidPointer(has, 4);
585e7bdbf8eSMatthew G. Knepley   *has = PETSC_FALSE;
586e7bdbf8eSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
5875cdeb108SMatthew G. Knepley   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
5885cdeb108SMatthew G. Knepley   if (exists) {
589e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
590*729ab6d9SBarry Smith     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
591e7bdbf8eSMatthew G. Knepley #else
592*729ab6d9SBarry Smith     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
593e7bdbf8eSMatthew G. Knepley #endif
594*729ab6d9SBarry Smith     if (dataset < 0) PetscFunctionReturn(0);
595*729ab6d9SBarry Smith     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
596*729ab6d9SBarry Smith     if (hhas < 0) {
597*729ab6d9SBarry Smith       PetscStackCallHDF5(H5Dclose,(dataset));
598*729ab6d9SBarry Smith       PetscFunctionReturn(0);
599*729ab6d9SBarry Smith     }
600*729ab6d9SBarry Smith     PetscStackCallHDF5(H5Dclose,(dataset));
6015cdeb108SMatthew G. Knepley     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
6025cdeb108SMatthew G. Knepley   }
603e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
604e7bdbf8eSMatthew G. Knepley }
605e7bdbf8eSMatthew G. Knepley 
606a75e6a4aSMatthew G. Knepley /*
607a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
608a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
609a75e6a4aSMatthew G. Knepley */
610a75e6a4aSMatthew G. Knepley static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
611a75e6a4aSMatthew G. Knepley 
612a75e6a4aSMatthew G. Knepley #undef __FUNCT__
613a75e6a4aSMatthew G. Knepley #define __FUNCT__ "PETSC_VIEWER_HDF5_"
614a75e6a4aSMatthew G. Knepley /*@C
615a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
616a75e6a4aSMatthew G. Knepley 
617a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
618a75e6a4aSMatthew G. Knepley 
619a75e6a4aSMatthew G. Knepley   Input Parameter:
620a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
621a75e6a4aSMatthew G. Knepley 
622a75e6a4aSMatthew G. Knepley   Level: intermediate
623a75e6a4aSMatthew G. Knepley 
624a75e6a4aSMatthew G. Knepley   Options Database Keys:
625a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
626a75e6a4aSMatthew G. Knepley 
627a75e6a4aSMatthew G. Knepley   Environmental variables:
628a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
629a75e6a4aSMatthew G. Knepley 
630a75e6a4aSMatthew G. Knepley   Notes:
631a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
632a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
633a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
634a75e6a4aSMatthew G. Knepley 
635a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
636a75e6a4aSMatthew G. Knepley @*/
637a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
638a75e6a4aSMatthew G. Knepley {
639a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
640a75e6a4aSMatthew G. Knepley   PetscBool      flg;
641a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
642a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
643a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
644a75e6a4aSMatthew G. Knepley 
645a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
646a75e6a4aSMatthew G. Knepley   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
647a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
648a75e6a4aSMatthew G. Knepley     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
649a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
650a75e6a4aSMatthew G. Knepley   }
651a75e6a4aSMatthew G. Knepley   ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
652a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
653a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
654a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
655a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
656a75e6a4aSMatthew G. Knepley     if (!flg) {
657a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
658a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
659a75e6a4aSMatthew G. Knepley     }
660a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
661a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
662a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
663a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
664a75e6a4aSMatthew G. Knepley     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
665a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
666a75e6a4aSMatthew G. Knepley   }
667a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
668a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
669a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
670a75e6a4aSMatthew G. Knepley }
671